CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/CalibTracker/SiPixelGainCalibration/plugins/SiPixelGainCalibrationAnalysis.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiPixelGainCalibrationAnalysis
00004 // Class:      SiPixelGainCalibrationAnalysis
00005 // 
00013 //
00014 // Original Author:  Freya Blekman
00015 //         Created:  Wed Nov 14 15:02:06 CET 2007
00016 // $Id: SiPixelGainCalibrationAnalysis.cc,v 1.41 2009/10/21 09:28:07 rougny Exp $
00017 //
00018 //
00019 
00020 // user include files
00021 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00022 
00023 #include "SiPixelGainCalibrationAnalysis.h"
00024 #include <sstream>
00025 #include <vector>
00026 #include <math.h>
00027 #include "TGraphErrors.h"
00028 #include "TMath.h"
00029 
00030 using std::cout;
00031 using std::endl;
00032 //
00033 // constructors and destructor
00034 //
00035 SiPixelGainCalibrationAnalysis::SiPixelGainCalibrationAnalysis(const edm::ParameterSet& iConfig):
00036   SiPixelOfflineCalibAnalysisBase(iConfig),
00037   conf_(iConfig),
00038   bookkeeper_(),
00039   bookkeeper_pixels_(),
00040   nfitparameters_(iConfig.getUntrackedParameter<int>("numberOfFitParameters",2)),
00041   fitfunction_(iConfig.getUntrackedParameter<std::string>("fitFunctionRootFormula","pol1")),
00042    listofdetids_(conf_.getUntrackedParameter<std::vector<uint32_t> >("listOfDetIDs")),
00043   ignoreMode_(conf_.getUntrackedParameter<bool>("ignoreMode",false)),
00044   reject_plateaupoints_(iConfig.getUntrackedParameter<bool>("suppressPlateauInFit",true)),
00045   reject_single_entries_(iConfig.getUntrackedParameter<bool>("suppressPointsWithOneEntryOrLess",true)),
00046   plateau_max_slope_(iConfig.getUntrackedParameter<double>("plateauSlopeMax",1.0)),
00047   reject_first_point_(iConfig.getUntrackedParameter<bool>("rejectVCalZero",true)),
00048   reject_badpoints_frac_(iConfig.getUntrackedParameter<double>("suppressZeroAndPlateausInFitFrac",0)),
00049   bookBIGCalibPayload_(iConfig.getUntrackedParameter<bool>("saveFullPayloads",false)),
00050   savePixelHists_(iConfig.getUntrackedParameter<bool>("savePixelLevelHists",false)),
00051   chi2Threshold_(iConfig.getUntrackedParameter<double>("minChi2NDFforHistSave",10)),
00052   chi2ProbThreshold_(iConfig.getUntrackedParameter<double>("minChi2ProbforHistSave",0.05)),
00053   maxGainInHist_(iConfig.getUntrackedParameter<double>("maxGainInHist",10)),
00054   maxChi2InHist_(iConfig.getUntrackedParameter<double>("maxChi2InHist",25)),
00055   saveALLHistograms_(iConfig.getUntrackedParameter<bool>("saveAllHistograms",false)),
00056  
00057 
00058   filldb_(iConfig.getUntrackedParameter<bool>("writeDatabase",false)),
00059   writeSummary_(iConfig.getUntrackedParameter<bool>("writeSummary",true)),
00060   recordName_(conf_.getParameter<std::string>("record")),
00061 
00062   appendMode_(conf_.getUntrackedParameter<bool>("appendMode",true)),
00063   /*theGainCalibrationDbInput_(0),
00064   theGainCalibrationDbInputOffline_(0),
00065   theGainCalibrationDbInputHLT_(0),
00066   theGainCalibrationDbInputService_(iConfig),*/
00067   gainlow_(10.),gainhi_(0.),pedlow_(255.),pedhi_(0.),
00068   useVcalHigh_(conf_.getParameter<bool>("useVCALHIGH")),
00069   scalarVcalHigh_VcalLow_(conf_.getParameter<double>("vcalHighToLowConversionFac"))
00070 {
00071   if(reject_single_entries_)
00072     min_nentries_=1;
00073   else
00074     min_nentries_=0;
00075   ::putenv((char*)"CORAL_AUTH_USER=me");
00076   ::putenv((char*)"CORAL_AUTH_PASSWORD=test");   
00077   edm::LogInfo("SiPixelGainCalibrationAnalysis") << "now using fit function " << fitfunction_ << ", which has " << nfitparameters_ << " free parameters. " << std::endl;
00078   func_= new TF1("func",fitfunction_.c_str(),0,256*scalarVcalHigh_VcalLow_);
00079   graph_ = new TGraphErrors();
00080   currentDetID_=0;
00081   summary_.open("SummaryPerDetID.txt");
00082   statusNumbers_ = new int[10];
00083   for(int ii=0;ii<10;ii++)
00084     statusNumbers_[ii]=0;
00085 }
00086 
00087 SiPixelGainCalibrationAnalysis::~SiPixelGainCalibrationAnalysis()
00088 {
00089 }
00090 // member functions
00091 //
00092 // ------------ method called once each job just before starting event loop  ------------
00093 
00094 std::vector<float> SiPixelGainCalibrationAnalysis::CalculateAveragePerColumn(uint32_t detid, std::string label){
00095   std::vector<float> result;
00096   int ncols= bookkeeper_[detid][label]->getNbinsX();
00097   int nrows= bookkeeper_[detid][label]->getNbinsY();
00098   for(int icol=1; icol<=ncols; ++icol){
00099     float val=0;
00100     float ntimes =0;
00101     for(int irow=1; irow<=nrows; ++irow){
00102       val+= bookkeeper_[detid][label]->getBinContent(icol,irow);
00103       ntimes++;
00104     }
00105     val/= ntimes;
00106     result.push_back(val);
00107   }
00108   return result;
00109 }
00110 
00111 bool
00112 SiPixelGainCalibrationAnalysis::checkCorrectCalibrationType()
00113 {
00114   if(calibrationMode_=="GainCalibration")
00115     return true;
00116   else if(ignoreMode_==true)
00117     return true;
00118   else if(calibrationMode_=="unknown"){
00119     edm::LogInfo("SiPixelGainCalibrationAnalysis") <<  "calibration mode is: " << calibrationMode_ << ", continuing anyway..." ;
00120     return true;
00121   }
00122   else{
00123     //    edm::LogError("SiPixelGainCalibrationAnalysis") << "unknown calibration mode for Gain calibration, should be \"Gain\" and is \"" << calibrationMode_ << "\"";
00124   }
00125   return false;
00126 }
00127 
00128 void SiPixelGainCalibrationAnalysis::calibrationSetup(const edm::EventSetup&)
00129 {
00130 }
00131 //------- summary printing method. Very verbose.
00132 void
00133 SiPixelGainCalibrationAnalysis::printSummary(){
00134 
00135   uint32_t detid=0;
00136   for(std::map<uint32_t,std::map<std::string,MonitorElement *> >::const_iterator idet = bookkeeper_.begin(); idet != bookkeeper_.end(); ++idet){
00137     if(detid==idet->first)
00138       continue;// only do things once per detid
00139     detid=idet->first;
00140     std::vector<float> gainvec=CalculateAveragePerColumn(detid,"gain_2d");
00141     std::vector<float> pedvec =CalculateAveragePerColumn(detid,"ped_2d");
00142     std::vector<float> chi2vec = CalculateAveragePerColumn(detid,"chi2_2d");
00143     std::ostringstream summarytext;
00144 
00145     summarytext << "Summary for det ID " << detid << "(" << translateDetIdToString(detid) << ")\n";
00146     summarytext << "\t Following: values per column: column #, gain, pedestal, chi2\n";
00147     for(uint32_t i=0; i<gainvec.size(); i++)
00148       summarytext << "\t " << i << " \t" << gainvec[i] << " \t" << pedvec[i] << " \t" << chi2vec[i] << "\n";
00149     summarytext << "\t list of pixels with high chi2 (chi2> " << chi2Threshold_ << "): \n";
00150 
00151     
00152     for(std::map<std::string, MonitorElement *>::const_iterator ipix = bookkeeper_pixels_[detid].begin(); ipix!=bookkeeper_pixels_[detid].end(); ++ipix)
00153       summarytext << "\t " << ipix->first << "\n";
00154     edm::LogInfo("SiPixelGainCalibrationAnalysis") << summarytext.str() << std::endl;
00155 
00156   }
00157   if(summary_.is_open()){
00158     summary_.close();
00159     summary_.open("Summary.txt");
00160     summary_<<"Total Number of Pixel computed :"<<statusNumbers_[9]<<endl;
00161     summary_<<"Number of pixel tagged with status :"<<endl;
00162     for(int ii=0;ii<9;ii++)
00163       summary_<<ii<<" -> "<<statusNumbers_[ii]<<" ~ "<<double(statusNumbers_[ii])/double(statusNumbers_[9])*100.<<" %"<<endl;
00164     
00165     summary_.close();
00166   
00167   }
00168 
00169 }
00170 
00171 // ------------ method called once each job just after ending the event loop  ------------
00172 
00173 void 
00174 SiPixelGainCalibrationAnalysis::calibrationEnd() {
00175 
00176   if(writeSummary_) printSummary();
00177   
00178   // this is where we loop over all histograms and save the database objects
00179   if(filldb_)
00180     fillDatabase();
00181 }
00182 //-----------method to fill the database
00183 void SiPixelGainCalibrationAnalysis::fillDatabase(){
00184  // only create when necessary.
00185   // process the minimum and maximum gain & ped values...
00186   edm::LogError("SiPixelGainCalibration::fillDatabase()") << "PLEASE do not fill the database directly from the gain calibration analyzer. This function is currently disabled and no DB payloads will be produced!" << std::endl;
00187 
00188 }
00189 // ------------ method called to do fits to all objects available  ------------
00190 bool
00191 SiPixelGainCalibrationAnalysis::doFits(uint32_t detid, std::vector<SiPixelCalibDigi>::const_iterator ipix)
00192 {
00193   float lowmeanval=255;
00194   float highmeanval=0;
00195   bool makehistopersistent = saveALLHistograms_;
00196   std::vector<uint32_t>::const_iterator detidfinder=find(listofdetids_.begin(),listofdetids_.end(),detid);
00197   if(detidfinder!=listofdetids_.end())
00198     makehistopersistent=true;
00199   // first, fill the input arrays to the TLinearFitter.
00200   double xvals[257];
00201   double yvals[256];
00202   double yerrvals[256];
00203   double xvalsall[257];
00204   double yvalsall[256];
00205   double yerrvalsall[256];
00206   int npoints=0;
00207   int nallpoints=0;
00208   bool use_point=true;
00209   int status=0;
00210   statusNumbers_[9]++;
00211   
00212   bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,0);
00213   if(writeSummary_ && detid!=currentDetID_){
00214     currentDetID_ = detid;
00215     summary_<<endl<<"DetId_"<<currentDetID_<<endl;
00216   }
00217   
00218   for(uint32_t ii=0; ii< ipix->getnpoints() && ii<200; ii++){
00219     //    std::cout << ipix->getsum(ii) << " " << ipix->getnentries(ii) << " " << ipix->getsumsquares(ii) << std::endl;
00220     nallpoints++;
00221     use_point=true;
00222     if(useVcalHigh_){
00223       xvalsall[ii]=vCalValues_[ii]*scalarVcalHigh_VcalLow_;
00224     }
00225     else
00226       xvalsall[ii]=vCalValues_[ii];
00227     yerrvalsall[ii]=yvalsall[ii]=0; 
00228   
00229     if(ipix->getnentries(ii)>min_nentries_){
00230       yvalsall[ii]=ipix->getsum(ii)/(float)ipix->getnentries(ii);
00231       yerrvalsall[ii]=ipix->getsumsquares(ii)/(float)(ipix->getnentries(ii));
00232       yerrvalsall[ii]-=pow(yvalsall[ii],2);
00233       yerrvalsall[ii]=sqrt(yerrvalsall[ii])/sqrt(ipix->getnentries(ii));
00234       
00235       if(yvalsall[ii]<lowmeanval)
00236         lowmeanval=yvalsall[ii];
00237       if(yvalsall[ii]>highmeanval)
00238         highmeanval=yvalsall[ii];
00239     }
00240   }
00241   
00242   // calculate plateau value from last 4 entries
00243   double plateauval=0;
00244   bool noPlateau=0;
00245   if(nallpoints>=4){
00246     for(int ii=nallpoints-1; ii>nallpoints-5; --ii) plateauval+=yvalsall[ii];
00247     plateauval/=4;
00248     for(int ii=nallpoints-1; ii>nallpoints-5; --ii){
00249       if(fabs(yvalsall[ii]-plateauval)>5){
00250         plateauval=255;
00251         noPlateau=1;
00252         continue;
00253       }
00254     }
00255     
00256     int NbofPointsInPlateau=0;
00257     for(int ii=0; ii<nallpoints; ++ii)
00258       if(fabs(yvalsall[ii]-plateauval)<10 || yvalsall[ii]==0) NbofPointsInPlateau++;
00259     //summary_<<"row_"<<ipix->row()<<" col_"<<ipix->col()<<"   "<<plateauval<<"  "<<NbofPointsInPlateau<<"  "<<nallpoints<<endl;
00260     if(NbofPointsInPlateau>=(nallpoints-2)){
00261       status=-2;
00262       bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,status);
00263       if(writeSummary_){
00264         summary_<<"row_"<<ipix->row()<<" col_"<<ipix->col()<<" status_"<<status<<endl;
00265         statusNumbers_[abs(status)]++;
00266       }
00267       return false;
00268     }
00269   }
00270   else plateauval=255;
00271     
00272   double maxgoodvalinfit=plateauval*(1.-reject_badpoints_frac_);
00273   npoints=0;
00274   for(int ii=0; ii<nallpoints; ++ii){
00275    
00276     // now selecting the appropriate points for the fit.
00277     use_point=true;
00278     if(reject_first_point_ && xvalsall[ii]<0.1)
00279       use_point=false;
00280     if(ipix->getnentries(ii)<=min_nentries_ && reject_single_entries_)
00281       use_point=false;
00282     if(ipix->getnentries(ii)==0 && reject_badpoints_)
00283       use_point=false;
00284     if(yvalsall[ii]>maxgoodvalinfit && !noPlateau)
00285       use_point=false;
00286     if(ii>1 && fabs(yvalsall[ii]-yvalsall[ii-1])<5. && yvalsall[ii]>0.8*maxgoodvalinfit && reject_plateaupoints_){
00287       use_point=false;
00288       break;
00289     }
00290     
00291     if(use_point){
00292       xvals[npoints]=xvalsall[ii];
00293       yvals[npoints]=yvalsall[ii];
00294       yerrvals[npoints]=yerrvalsall[ii];
00295       npoints++;
00296     }
00297   }
00298   
00299   float chi2,slope,intercept,prob,slopeerror,intercepterror;
00300   prob=chi2=-1;
00301   slope=intercept=slopeerror=intercepterror=0;
00302   
00303 
00304   // now check on number of points. If bad just start taking the first 4:
00305 
00306   if(npoints<4){
00307     npoints=0;
00308     for(int ii=0; ii<nallpoints && npoints<4 && yvalsall[ii]<plateauval*0.97; ++ii){
00309       if(yvalsall[ii]>0){
00310         if(ii>0 && yvalsall[ii]-yvalsall[ii-1]<0.1)
00311           continue;
00312         xvals[npoints]=xvalsall[ii];
00313         yvals[npoints]=yvalsall[ii];
00314         yerrvals[npoints]=yerrvalsall[ii];
00315         npoints++;
00316       }
00317     }
00318   }
00319   if(npoints<2){
00320     status = -7;
00321     bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,status);
00322     if(writeSummary_){
00323       summary_<<"row_"<<ipix->row()<<" col_"<<ipix->col()<<" status_"<<status<<endl;
00324       statusNumbers_[abs(status)]++;
00325     }
00326     std::ostringstream pixelinfo;
00327     pixelinfo << "GainCurve_row_" << ipix->row() << "_col_" << ipix->col();
00328     std::string tempname=translateDetIdToString(detid);
00329     tempname+="_";
00330     tempname+=pixelinfo.str();
00331     setDQMDirectory(detid);
00332     bookkeeper_pixels_[detid][pixelinfo.str()] = bookDQMHistogram1D(detid,pixelinfo.str(),tempname,105*nallpoints,xvalsall[0],xvalsall[nallpoints-1]*1.05);
00333     for(int ii=0; ii<nallpoints; ++ii)
00334       bookkeeper_pixels_[detid][pixelinfo.str()]->Fill(xvalsall[ii],yvalsall[ii]);
00335     return false;
00336   }
00337     
00338   //  std::cout << "starting fit!" << std::endl;
00339   graph_->Set(npoints);
00340   
00341   func_->SetParameter(0,50.);
00342   func_->SetParameter(1,0.25);
00343   for(int ipointtemp=0; ipointtemp<npoints; ++ipointtemp){
00344     graph_->SetPoint(ipointtemp,xvals[ipointtemp],yvals[ipointtemp]);
00345     graph_->SetPointError(ipointtemp,0,yerrvals[ipointtemp]);
00346   }
00347   Int_t tempresult = graph_->Fit("func","FQ0N");
00348   slope=func_->GetParameter(1);
00349   slopeerror=func_->GetParError(1);
00350   intercept=func_->GetParameter(0);
00351   intercepterror=func_->GetParError(0);
00352   chi2=func_->GetChisquare()/((float)npoints-func_->GetNpar());
00353   prob= TMath::Prob(func_->GetChisquare(),npoints-func_->GetNpar());
00354   size_t ntimes=0;
00355   while((isnan(slope) || isnan(intercept) )&& ntimes<10){
00356     ntimes++;
00357     makehistopersistent=true;
00358     //    std::cout << slope << " " << intercept << " " << prob << std::endl;
00359     edm::LogWarning("SiPixelGainCalibrationAnalysis") << "impossible to fit values, try " << ntimes << ": " ;
00360     for(int ii=0; ii<npoints; ++ii){
00361       edm::LogWarning("SiPixelGainCalibrationAnalysis")<< "vcal " << xvals[ii] << " response: " << yvals[ii] << "+/-" << yerrvals[ii] << std::endl; 
00362     } 
00363     tempresult = graph_->Fit("func","FQ0NW");
00364     slope=func_->GetParameter(1);
00365     slopeerror=func_->GetParError(1);
00366     intercept = func_->GetParameter(0);
00367     intercepterror=func_->GetParError(0);
00368     chi2=func_->GetChisquare()/((float)npoints-func_->GetNpar());
00369     prob= TMath::Prob(func_->GetChisquare(),npoints-func_->GetNpar());
00370   }
00371   
00372   if(tempresult==0)
00373     status=1;
00374   else
00375     status=0;
00376   if(slope!=0)
00377     slope = 1./slope;
00378   if(isnan(slope) || isnan(intercept)){
00379     status=-6;
00380     bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,status);
00381     if(writeSummary_){
00382       summary_<<"row_"<<ipix->row()<<" col_"<<ipix->col()<<" status_"<<status<<endl;
00383       statusNumbers_[abs(status)]++;
00384     }
00385     //return false;
00386   }
00387   if(chi2>chi2Threshold_ && chi2Threshold_>=0)
00388     status=5;
00389   if(prob<chi2ProbThreshold_)
00390     status=5;
00391   if(noPlateau)
00392     status=3;
00393   if(nallpoints<4)
00394     status=-7;
00395   if(TMath::Abs(slope>maxGainInHist_) || slope < 0)
00396     status=-8;
00397   if(status!=1)
00398     makehistopersistent=true;
00399   statusNumbers_[abs(status)]++;
00400 
00401 
00402   if(slope<gainlow_)
00403     gainlow_=slope;
00404   if(slope>gainhi_)
00405     gainhi_=slope;
00406   if(intercept>pedhi_)
00407     pedhi_=intercept;
00408   if(intercept<pedlow_)
00409     pedlow_=intercept;
00410   bookkeeper_[detid]["gain_1d"]->Fill(slope);
00411   if(slope>maxGainInHist_){
00412     makehistopersistent=true;
00413     edm::LogWarning("SiPixelGainCalibration") << "For DETID " << detid << "pixel row,col " << ipix->row() << "," << ipix->col() << " Gain was measured to be " << slope << " which is outside the range of the summary plot (" <<maxGainInHist_ << ") !!!! " << std::endl;
00414   }
00415   bookkeeper_[detid]["dynamicrange_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,highmeanval-lowmeanval);
00416   bookkeeper_[detid]["plateau_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,highmeanval);
00417   bookkeeper_[detid]["gain_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,slope);
00418   bookkeeper_[detid]["errorgain_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,slopeerror);
00419   bookkeeper_[detid]["ped_1d"]->Fill(intercept);
00420   bookkeeper_[detid]["ped_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,intercept);
00421   bookkeeper_[detid]["errorped_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,intercepterror);
00422   bookkeeper_[detid]["chi2_1d"]->Fill(chi2);
00423   bookkeeper_[detid]["chi2_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,chi2);
00424   bookkeeper_[detid]["prob_1d"]->Fill(prob);
00425   bookkeeper_[detid]["prob_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,prob);
00426   bookkeeper_[detid]["lowpoint_1d"]->Fill(xvals[0]);
00427   bookkeeper_[detid]["lowpoint_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,xvals[0]);
00428   bookkeeper_[detid]["highpoint_1d"]->Fill(xvals[npoints-1]);
00429   bookkeeper_[detid]["highpoint_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,xvals[npoints-1]);
00430   bookkeeper_[detid]["nfitpoints_1d"]->Fill(npoints);
00431   bookkeeper_[detid]["endpoint_1d"]->Fill((255 - intercept)*slope);
00432   bookkeeper_[detid]["status_2d"]->setBinContent(ipix->col()+1,ipix->row()+1,status);
00433 
00434     if(!savePixelHists_)
00435     return true;
00436   if(detidfinder==listofdetids_.end() && listofdetids_.size()!=0)
00437     return true;
00438   if(makehistopersistent){
00439     std::ostringstream pixelinfo;
00440     pixelinfo << "GainCurve_row_" << ipix->row() << "_col_" << ipix->col();
00441     std::string tempname=translateDetIdToString(detid);
00442     tempname+="_";
00443     tempname+=pixelinfo.str();
00444     
00445     // and book the histo
00446     // fill the last value of the vcal array...   
00447 
00448     setDQMDirectory(detid);
00449     bookkeeper_pixels_[detid][pixelinfo.str()] =  bookDQMHistogram1D(detid,pixelinfo.str(),tempname,105*nallpoints,xvalsall[0],xvalsall[nallpoints-1]*1.05);
00450     
00451     edm::LogInfo("SiPixelGainCalibrationAnalysis") << "now saving histogram for pixel " << tempname << ", gain = " << slope << ", pedestal = " << intercept << ", chi2/NDF=" << chi2 << "(prob:" << prob << "), fit status " << status;
00452     for(int ii=0; ii<nallpoints; ++ii){
00453       //      std::cout << xvalsall[ii]<<","<<yvalsall[ii]<< " " << tempfloats[ii] << std::endl;
00454       bookkeeper_pixels_[detid][pixelinfo.str()]->Fill(xvalsall[ii],yvalsall[ii]);
00455     }
00456     
00457     //    addTF1ToDQMMonitoringElement(bookkeeper_pixels_[detid][pixelinfo.str()],func_);
00458     
00459     
00460     if(writeSummary_){
00461       summary_<<"row_"<<ipix->row()<<" col_"<<ipix->col();
00462       summary_<<" status_"<<status;
00463       summary_<<endl;
00464       
00465       //std::cout<<detid<<"  "<<"row " <<ipix->row()<<" col "<<ipix->col()<<"  "<<status<<"  "<<chi2<<"  "<<prob<<"  "<<npoints<<"  "<<xvals[0]<<"  "<<xvals[npoints-1]<<"  "<<plateauval<<std::endl;
00466     }
00467   } 
00468   return true;
00469 }
00470 // ------------ method called to do fill new detids  ------------
00471 void 
00472 SiPixelGainCalibrationAnalysis::newDetID(uint32_t detid)
00473 {
00474   setDQMDirectory(detid);
00475   std::string tempname=translateDetIdToString(detid);
00476   bookkeeper_[detid]["gain_1d"] = bookDQMHistogram1D(detid,"Gain1d","gain for "+tempname,100,0.,maxGainInHist_);
00477   bookkeeper_[detid]["gain_2d"] = bookDQMHistoPlaquetteSummary2D(detid, "Gain2d","gain for "+tempname);
00478   bookkeeper_[detid]["errorgain_2d"] = bookDQMHistoPlaquetteSummary2D(detid, "ErrorGain2d","error on gain for "+tempname);
00479   bookkeeper_[detid]["ped_1d"] = bookDQMHistogram1D(detid,"Pedestal1d","pedestal for "+tempname,256,0.,256.0);
00480   bookkeeper_[detid]["ped_2d"] = bookDQMHistoPlaquetteSummary2D(detid,"Pedestal2d","pedestal for "+tempname);
00481   bookkeeper_[detid]["errorped_2d"] = bookDQMHistoPlaquetteSummary2D(detid,"ErrorPedestal2d","error on pedestal for "+tempname);
00482   bookkeeper_[detid]["chi2_1d"] = bookDQMHistogram1D(detid,"GainChi2NDF1d","#chi^{2}/NDOF for "+tempname,100,0.,maxChi2InHist_);
00483   bookkeeper_[detid]["chi2_2d"] = bookDQMHistoPlaquetteSummary2D(detid,"GainChi2NDF2d","#chi^{2}/NDOF for "+tempname);
00484   bookkeeper_[detid]["prob_1d"] = bookDQMHistogram1D(detid,"GainChi2Prob1d","P(#chi^{2},NDOF) for "+tempname,100,0.,1.0);
00485   bookkeeper_[detid]["prob_2d"] = bookDQMHistoPlaquetteSummary2D(detid,"GainChi2Prob2d","P(#chi^{2},NDOF) for "+tempname);
00486   bookkeeper_[detid]["status_2d"] = bookDQMHistoPlaquetteSummary2D(detid,"GainFitResult2d","Fit result for "+tempname);
00487   bookkeeper_[detid]["endpoint_1d"]= bookDQMHistogram1D(detid,"GainEndPoint1d","point where fit meets ADC=255 for "+tempname,256,0.,256.*scalarVcalHigh_VcalLow_);
00488   bookkeeper_[detid]["lowpoint_1d"]= bookDQMHistogram1D(detid,"GainLowPoint1d","lowest fit point for "+tempname,256,0.,256.*scalarVcalHigh_VcalLow_);
00489   bookkeeper_[detid]["highpoint_1d"]= bookDQMHistogram1D(detid,"GainHighPoint1d","highest fit point for "+tempname,256,0.,256.*scalarVcalHigh_VcalLow_);
00490   bookkeeper_[detid]["nfitpoints_1d"]= bookDQMHistogram1D(detid,"GainNPoints1d","number of fit point for "+tempname,20,0.,20);
00491   bookkeeper_[detid]["dynamicrange_2d"]=bookDQMHistoPlaquetteSummary2D(detid,"GainDynamicRange2d","Difference lowest and highest points on gain curve for "+tempname); 
00492   bookkeeper_[detid]["lowpoint_2d"]=bookDQMHistoPlaquetteSummary2D(detid,"GainLowPoint2d","lowest fit point for "+tempname);
00493   bookkeeper_[detid]["highpoint_2d"]=bookDQMHistoPlaquetteSummary2D(detid,"GainHighPoint2d","highest fit point for "+tempname);
00494   bookkeeper_[detid]["plateau_2d"]=bookDQMHistoPlaquetteSummary2D(detid,"GainSaturate2d","Highest points on gain curve for "+tempname);
00495 
00496 }
00497 //define this as a plug-in
00498 DEFINE_FWK_MODULE(SiPixelGainCalibrationAnalysis);