Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:03:24

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiPixelGainCalibrationReadDQMFile
0004 // Class:      SiPixelGainCalibrationReadDQMFile
0005 //
0006 /**\class SiPixelGainCalibrationReadDQMFile SiPixelGainCalibrationReadDQMFile.cc CalibTracker/SiPixelGainCalibrationReadDQMFile/src/SiPixelGainCalibrationReadDQMFile.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Freya BLEKMAN
0015 //         Created:  Tue Aug  5 16:22:46 CEST 2008
0016 // $Id: SiPixelGainCalibrationReadDQMFile.cc,v 1.7 2010/01/12 11:29:54 rougny Exp $
0017 //
0018 //
0019 
0020 // system include files
0021 #include <memory>
0022 #include <fstream>
0023 #include <sys/stat.h>
0024 
0025 // user include files
0026 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationService.h"
0027 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0028 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0029 #include "CondFormats/SiPixelObjects/interface/SiPixelCalibConfiguration.h"
0030 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibration.h"
0031 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationForHLT.h"
0032 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationOffline.h"
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include "FWCore/Framework/interface/Frameworkfwd.h"
0035 #include "FWCore/Framework/interface/MakerMacros.h"
0036 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0038 #include "FWCore/ServiceRegistry/interface/Service.h"
0039 #include "FWCore/ServiceRegistry/interface/Service.h"
0040 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0041 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0042 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0043 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0044 
0045 // ROOT includes
0046 #include "TDirectory.h"
0047 #include "TFile.h"
0048 #include "TH2F.h"
0049 #include "TKey.h"
0050 #include "TList.h"
0051 #include "TString.h"
0052 #include "TTree.h"
0053 
0054 //
0055 // class decleration
0056 //
0057 
0058 class SiPixelGainCalibrationReadDQMFile : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0059 public:
0060   explicit SiPixelGainCalibrationReadDQMFile(const edm::ParameterSet &);
0061 
0062 private:
0063   void analyze(const edm::Event &, const edm::EventSetup &) final;
0064   // functions added by F.B.
0065   void fillDatabase(const edm::EventSetup &iSetup, TFile *);
0066   std::unique_ptr<TFile> getHistograms();
0067   // ----------member data ---------------------------
0068   std::map<uint32_t, std::map<std::string, TString> > bookkeeper_;
0069   std::map<uint32_t, std::map<double, double> > Meankeeper_;
0070   std::map<uint32_t, std::vector<std::map<int, int> > > noisyPixelsKeeper_;
0071 
0072   const bool appendMode_;
0073   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> pddToken_;
0074   SiPixelGainCalibrationService theGainCalibrationDbInputService_;
0075   std::unique_ptr<TH2F> defaultGain_;
0076   std::unique_ptr<TH2F> defaultPed_;
0077   std::unique_ptr<TH2F> defaultChi2_;
0078   std::unique_ptr<TH2F> defaultFitResult_;
0079   std::unique_ptr<TH1F> meanGainHist_;
0080   std::unique_ptr<TH1F> meanPedHist_;
0081   const std::string record_;
0082   // keep track of lowest and highest vals for range
0083   float gainlow_;
0084   float gainhi_;
0085   float pedlow_;
0086   float pedhi_;
0087   const bool usemeanwhenempty_;
0088   const std::string rootfilestring_;
0089   float gainmax_;
0090   float pedmax_;
0091   const double badchi2_;
0092   const size_t nmaxcols;
0093   const size_t nmaxrows;
0094 };
0095 
0096 void SiPixelGainCalibrationReadDQMFile::fillDatabase(const edm::EventSetup &iSetup, TFile *therootfile) {
0097   // only create when necessary.
0098   // process the minimum and maximum gain & ped values...
0099   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now starting db fill!!!" << std::endl;
0100 
0101   std::map<uint32_t, std::pair<TString, int> > badresults;
0102 
0103   edm::Service<TFileService> fs;
0104 
0105   //   bool usedmean=false;
0106   //   bool isdead=false;
0107   TH1F *VCAL_endpoint =
0108       fs->make<TH1F>("VCAL_endpoint", "value where response = 255 ( x = (255 - ped)/gain )", 256, 0, 256);
0109   TH1F *goodgains = fs->make<TH1F>("goodgains", "gain values", 100, 0, 10);
0110   TH1F *goodpeds = fs->make<TH1F>("goodpeds", "pedestal values", 356, -100, 256);
0111   TH1F *totgains = fs->make<TH1F>("totgains", "gain values", 200, 0, 10);
0112   TH1F *totpeds = fs->make<TH1F>("totpeds", "pedestal values", 356, -100, 256);
0113   TTree *tree = new TTree("tree", "tree");
0114   int detidfortree, rowfortree, colfortree, useddefaultfortree;
0115   float pedfortree, gainfortree, chi2fortree;
0116   tree->Branch("detid", &detidfortree, "detid/I");
0117   tree->Branch("row", &rowfortree, "row/I");
0118   tree->Branch("col", &colfortree, "col/I");
0119   tree->Branch("defaultval", &useddefaultfortree, "defaultval/I");
0120   tree->Branch("ped", &pedfortree, "ped/F");
0121   tree->Branch("gain", &gainfortree, "gain/F");
0122   tree->Branch("chi2", &chi2fortree, "chi2/F");
0123 
0124   //   TH1F *gainPerDetid;
0125   //   TH1F *pedPerDetid;
0126 
0127   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now filling record " << record_ << std::endl;
0128   if (record_ != "SiPixelGainCalibrationForHLTRcd" && record_ != "SiPixelGainCalibrationOfflineRcd") {
0129     edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0130         << "you passed record " << record_ << ", which I have no idea what to do with!" << std::endl;
0131     return;
0132   }
0133   if (gainlow_ > gainhi_) {
0134     float temp = gainhi_;
0135     gainhi_ = gainlow_;
0136     gainlow_ = temp;
0137   }
0138   if (pedlow_ > pedhi_) {
0139     float temp = pedhi_;
0140     pedhi_ = pedlow_;
0141     pedlow_ = temp;
0142   }
0143   if (gainhi_ > gainmax_)
0144     gainhi_ = gainmax_;
0145   if (pedhi_ > pedmax_)
0146     pedhi_ = pedmax_;
0147   float badpedval = pedlow_ - 200;
0148   float badgainval = gainlow_ - 200;
0149   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now filling db: values: pedlow, hi: " << pedlow_ << ", "
0150                                                      << pedhi_ << " and gainlow, hi: " << gainlow_ << ", " << gainhi_;
0151   float meangain = meanGainHist_->GetMean();
0152   float meanped = meanPedHist_->GetMean();
0153   //  edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << ", and mean gain " << meangain<< ", ped " << meanped ;
0154   //  edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << std::endl;
0155 
0156   // and fill the dummy histos:
0157 
0158   for (size_t icol = 0; icol < nmaxcols; ++icol) {
0159     for (size_t irow = 0; irow < nmaxrows; ++irow) {
0160       defaultGain_->SetBinContent(icol + 1, irow + 1, meangain);
0161       defaultPed_->SetBinContent(icol + 1, irow + 1, meanped);
0162       defaultChi2_->SetBinContent(icol + 1, irow + 1, 1.0);
0163       defaultFitResult_->SetBinContent(icol + 1, irow + 1, 0);
0164     }
0165   }
0166 
0167   SiPixelGainCalibration theGainCalibrationDbInput(pedlow_ * 0.999, pedhi_ * 1.001, gainlow_ * 0.999, gainhi_ * 1.001);
0168   SiPixelGainCalibrationForHLT theGainCalibrationDbInputHLT(
0169       pedlow_ * 0.999, pedhi_ * 1.001, gainlow_ * 0.999, gainhi_ * 1.001);
0170   SiPixelGainCalibrationOffline theGainCalibrationDbInputOffline(
0171       pedlow_ * 0.999, pedhi_ * 1.001, gainlow_ * 0.999, gainhi_ * 1.001);
0172 
0173   uint32_t nchannels = 0;
0174   uint32_t nmodules = 0;
0175   edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0176       << "now starting loop on detids, there are " << bookkeeper_.size() << " histograms to consider..." << std::endl;
0177   uint32_t detid = 0;
0178   therootfile->cd();
0179   const TrackerGeometry *pDD = &iSetup.getData(pddToken_);
0180   edm::LogInfo("SiPixelCondObjOfflineBuilder") << " There are " << pDD->dets().size() << " detectors" << std::endl;
0181 
0182   for (TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++) {
0183     detid = 0;
0184     if (dynamic_cast<PixelGeomDetUnit const *>((*it)) != nullptr)
0185       detid = ((*it)->geographicalId()).rawId();
0186     if (detid == 0)
0187       continue;
0188     //if(detid!=344076812) continue;
0189     int badDetId = 0;
0190     //if(detid==302123296||detid==302126596) badDetId=1;;
0191     //if(detid!=302058516) continue;
0192     useddefaultfortree = 0;
0193     edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0194         << "now creating database object for detid " << detid
0195         << std::endl;  //<< " " << bookkeeper_[detid]["gain_2d"] << " " << bookkeeper_[detid]["ped_2d"] << std::endl; //edm::LogPrint("SiPixelGainCalibrationReadDQMFile")<< " nrows:" << nrows << " ncols: " << ncols << std::endl;
0196 
0197     // Get the module sizes.
0198     TH2F *tempchi2;
0199     TH2F *tempfitresult;
0200     TH2F *tempgain;
0201     TH2F *tempped;
0202     TString tempgainstring;
0203 
0204     if (!badDetId) {
0205       TString tempchi2string = bookkeeper_[detid]["chi2prob_2d"];
0206       tempchi2 = dynamic_cast<TH2F *>(therootfile->Get(tempchi2string));
0207       if (tempchi2 == nullptr || badDetId) {
0208         tempchi2 = defaultChi2_.get();
0209         useddefaultfortree = 1;
0210       }
0211       TString tempfitresultstring = bookkeeper_[detid]["fitresult_2d"];
0212       tempfitresult = dynamic_cast<TH2F *>(therootfile->Get(tempfitresultstring));
0213       if (tempfitresult == nullptr) {
0214         tempfitresult = defaultFitResult_.get();
0215         useddefaultfortree = 1;
0216       }
0217       TString tempgainstring = bookkeeper_[detid]["gain_2d"];
0218       tempgain = dynamic_cast<TH2F *>(therootfile->Get(tempgainstring));
0219       if (tempgain == nullptr) {
0220         //  edm::LogPrint("SiPixelGainCalibrationReadDQMFile") <<"WARNING, gain histo " << bookkeeper_[detid]["gain_2d"] << " does not exist, using default instead" << std::endl;
0221         tempgain = defaultGain_.get();
0222         useddefaultfortree = 1;
0223       }
0224       TString temppedstring = bookkeeper_[detid]["ped_2d"];
0225       tempped = dynamic_cast<TH2F *>(therootfile->Get(temppedstring));
0226       if (tempped == nullptr) {
0227         //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") <<"WARNING, ped histo " << bookkeeper_[detid]["ped_2d"] << " for detid " << detid << " does not exist, using default instead" << std::endl;
0228         std::pair<TString, int> tempval(tempgainstring, 0);
0229         badresults[detid] = tempval;
0230         tempped = defaultPed_.get();
0231         useddefaultfortree = 1;
0232       }
0233     } else {
0234       tempchi2 = defaultChi2_.get();
0235       tempgain = defaultGain_.get();
0236       tempfitresult = defaultFitResult_.get();
0237       std::pair<TString, int> tempval(tempgainstring, 0);
0238       badresults[detid] = tempval;
0239       tempped = defaultPed_.get();
0240       useddefaultfortree = 1;
0241     }
0242 
0243     const PixelGeomDetUnit *pixDet = dynamic_cast<const PixelGeomDetUnit *>((*it));
0244     const PixelTopology &topol = pixDet->specificTopology();
0245     // Get the module sizes.
0246     size_t nrows = topol.nrows();     // rows in x
0247     size_t ncols = topol.ncolumns();  // cols in y
0248 
0249     //    int nrows=tempgain->GetNbinsY();
0250     //    int ncols=tempgain->GetNbinsX();
0251     //    edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "next histo " << tempgain->GetTitle() << " has nrow,ncol:" << nrows << ","<< ncols << std::endl;
0252     size_t nrowsrocsplit = theGainCalibrationDbInputHLT.getNumberOfRowsToAverageOver();
0253     if (theGainCalibrationDbInputOffline.getNumberOfRowsToAverageOver() != nrowsrocsplit)
0254       throw cms::Exception("GainCalibration Payload configuration error")
0255           << "[SiPixelGainCalibrationAnalysis::fillDatabase] ERROR the SiPixelGainCalibrationOffline and "
0256              "SiPixelGainCalibrationForHLT database payloads have different settings for the number of rows per roc: "
0257           << theGainCalibrationDbInputHLT.getNumberOfRowsToAverageOver() << "(HLT), "
0258           << theGainCalibrationDbInputOffline.getNumberOfRowsToAverageOver() << "(offline)";
0259     std::vector<char> theSiPixelGainCalibrationPerPixel;
0260     std::vector<char> theSiPixelGainCalibrationPerColumn;
0261     std::vector<char> theSiPixelGainCalibrationGainPerColPedPerPixel;
0262 
0263     //Get mean of gain/pedestal of this Detid
0264     meangain = meanGainHist_->GetMean();
0265     meanped = meanPedHist_->GetMean();
0266     int npix = 0;
0267     double meanGainForThisModule = 0;
0268     double meanPedForThisModule = 0;
0269     for (size_t icol = 1; icol <= ncols; icol++) {
0270       for (size_t jrow = 1; jrow <= nrows; jrow++) {
0271         if (tempfitresult->GetBinContent(icol, jrow) > 0) {
0272           npix++;
0273           meanGainForThisModule += tempgain->GetBinContent(icol, jrow);
0274           meanPedForThisModule += tempped->GetBinContent(icol, jrow);
0275         }
0276       }
0277     }
0278     if (npix != 0)
0279       meanPedForThisModule /= npix;
0280     if (npix != 0)
0281       meanGainForThisModule /= npix;
0282     if (usemeanwhenempty_) {
0283       if (meanGainForThisModule > gainlow_ && meanGainForThisModule < gainhi_ && npix > 100)
0284         meangain = meanGainForThisModule;
0285       if (meanPedForThisModule > pedlow_ && meanPedForThisModule < pedhi_ && npix > 100)
0286         meanped = meanPedForThisModule;
0287     }
0288 
0289     // Loop over columns and rows of this DetID
0290     float peds[160];
0291     float gains[160];
0292     float pedforthiscol[2];
0293     float gainforthiscol[2];
0294     int nusedrows[2];
0295     //edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << pedlow_<<"  "<<pedhi_<<"  "<<gainlow_<<"  "<<gainhi_<< std::endl;
0296     for (size_t icol = 1; icol <= ncols; icol++) {
0297       nusedrows[0] = nusedrows[1] = 0;
0298       pedforthiscol[0] = pedforthiscol[1] = 0;
0299       gainforthiscol[0] = gainforthiscol[1] = 0;
0300       //   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now lookign at col " << icol << std::endl;
0301       for (size_t jrow = 1; jrow <= nrows; jrow++) {
0302         size_t iglobalrow = 0;
0303         if (jrow > nrowsrocsplit)
0304           iglobalrow = 1;
0305         peds[jrow] = badpedval;
0306         gains[jrow] = badgainval;
0307         float ped = tempped->GetBinContent(icol, jrow);
0308         float gain = tempgain->GetBinContent(icol, jrow);
0309         float chi2 = tempchi2->GetBinContent(icol, jrow);
0310         float fitresult = tempfitresult->GetBinContent(icol, jrow);
0311 
0312         if (ped > pedlow_ && gain > gainlow_ && ped < pedhi_ && gain < gainhi_ && (fitresult > 0)) {
0313           VCAL_endpoint->Fill((255 - ped) / gain);
0314           peds[jrow] = ped;
0315           gains[jrow] = gain;
0316           pedforthiscol[iglobalrow] += ped;
0317           gainforthiscol[iglobalrow] += gain;
0318           nusedrows[iglobalrow]++;
0319           goodpeds->Fill(ped);
0320           goodgains->Fill(gain);
0321           detidfortree = detid;
0322           rowfortree = jrow - 1;
0323           colfortree = icol - 1;
0324           gainfortree = gain;
0325           pedfortree = ped;
0326           chi2fortree = chi2;
0327           tree->Fill();
0328         } else {
0329           if (usemeanwhenempty_) {
0330             peds[jrow] = meanped;
0331             gains[jrow] = meangain;
0332             std::pair<TString, int> tempval(tempgainstring, 1);
0333             badresults[detid] = tempval;
0334           } else {
0335             std::pair<TString, int> tempval(tempgainstring, 2);
0336             badresults[detid] = tempval;
0337             // if everything else fails: set the gain & ped now to dead
0338             peds[jrow] = badpedval;
0339             gains[jrow] = badgainval;
0340           }
0341         }
0342         //edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << detid << " " << jrow << " " << icol << " " << peds[jrow] << " " << gains[jrow] << " " << chi2 << " " << fitresult << std::endl;
0343 
0344         totgains->Fill(gains[jrow]);
0345         totpeds->Fill(peds[jrow]);
0346 
0347       }  // now collected all info, actually do the filling
0348 
0349       for (size_t jrow = 1; jrow <= nrows; jrow++) {
0350         nchannels++;
0351         size_t iglobalrow = 0;
0352         if (jrow > nrowsrocsplit)
0353           iglobalrow = 1;
0354         float ped = peds[jrow];
0355         float gain = gains[jrow];
0356 
0357         if (ped > pedlow_ && gain > gainlow_ && ped < pedhi_ && gain < gainhi_) {
0358           theGainCalibrationDbInput.setData(ped, gain, theSiPixelGainCalibrationPerPixel);
0359           theGainCalibrationDbInputOffline.setDataPedestal(ped, theSiPixelGainCalibrationGainPerColPedPerPixel);
0360         } else {
0361           theGainCalibrationDbInput.setDeadPixel(theSiPixelGainCalibrationPerPixel);
0362           theGainCalibrationDbInputOffline.setDeadPixel(theSiPixelGainCalibrationGainPerColPedPerPixel);
0363         }
0364 
0365         if (jrow % nrowsrocsplit == 0) {
0366           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now in col " << icol << " " << jrow << " " << iglobalrow << std::endl;
0367           if (nusedrows[iglobalrow] > 0) {
0368             pedforthiscol[iglobalrow] /= (float)nusedrows[iglobalrow];
0369             gainforthiscol[iglobalrow] /= (float)nusedrows[iglobalrow];
0370           }
0371           if (gainforthiscol[iglobalrow] > gainlow_ && gainforthiscol[iglobalrow] < gainhi_ &&
0372               pedforthiscol[iglobalrow] > pedlow_ && pedforthiscol[iglobalrow] < pedhi_) {  // good
0373             //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "setting ped & col aves: " << pedforthiscol[iglobalrow] << " " <<  gainforthiscol[iglobalrow]<< std::endl;
0374           } else {
0375             if (usemeanwhenempty_) {
0376               //          edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "setting ped & col aves: " << pedforthiscol[iglobalrow] << " " <<  gainforthiscol[iglobalrow]<< std::endl;
0377               pedforthiscol[iglobalrow] = meanped;
0378               gainforthiscol[iglobalrow] = meangain;
0379               std::pair<TString, int> tempval(tempgainstring, 3);
0380               badresults[detid] = tempval;
0381             } else {  //make dead
0382               pedforthiscol[iglobalrow] = badpedval;
0383               gainforthiscol[iglobalrow] = badgainval;
0384               std::pair<TString, int> tempval(tempgainstring, 4);
0385               badresults[detid] = tempval;
0386             }
0387           }
0388 
0389           if (gainforthiscol[iglobalrow] > gainlow_ && gainforthiscol[iglobalrow] < gainhi_ &&
0390               pedforthiscol[iglobalrow] > pedlow_ && pedforthiscol[iglobalrow] < pedhi_) {
0391             theGainCalibrationDbInputOffline.setDataGain(
0392                 gainforthiscol[iglobalrow], nrowsrocsplit, theSiPixelGainCalibrationGainPerColPedPerPixel);
0393             theGainCalibrationDbInputHLT.setData(
0394                 pedforthiscol[iglobalrow], gainforthiscol[iglobalrow], theSiPixelGainCalibrationPerColumn);
0395           } else {
0396             //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << pedforthiscol[iglobalrow] << " " << gainforthiscol[iglobalrow] << std::endl;
0397             theGainCalibrationDbInputOffline.setDeadColumn(nrowsrocsplit,
0398                                                            theSiPixelGainCalibrationGainPerColPedPerPixel);
0399             theGainCalibrationDbInputHLT.setDeadColumn(nrowsrocsplit, theSiPixelGainCalibrationPerColumn);
0400           }
0401         }
0402       }
0403     }
0404 
0405     //    edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "setting range..." << std::endl;
0406     SiPixelGainCalibration::Range range(theSiPixelGainCalibrationPerPixel.begin(),
0407                                         theSiPixelGainCalibrationPerPixel.end());
0408     SiPixelGainCalibrationForHLT::Range hltrange(theSiPixelGainCalibrationPerColumn.begin(),
0409                                                  theSiPixelGainCalibrationPerColumn.end());
0410     SiPixelGainCalibrationOffline::Range offlinerange(theSiPixelGainCalibrationGainPerColPedPerPixel.begin(),
0411                                                       theSiPixelGainCalibrationGainPerColPedPerPixel.end());
0412 
0413     //    edm::LogPrint("SiPixelGainCalibrationReadDQMFile") <<"putting things in db..." << std::endl;
0414     // now start creating the various database objects
0415     if (!theGainCalibrationDbInput.put(detid, range, ncols))
0416       edm::LogError("SiPixelGainCalibrationAnalysis")
0417           << "warning: detid already exists for Offline (gain per col, ped per pixel) calibration database"
0418           << std::endl;
0419     if (!theGainCalibrationDbInputOffline.put(detid, offlinerange, ncols))
0420       edm::LogError("SiPixelGainCalibrationAnalysis")
0421           << "warning: detid already exists for Offline (gain per col, ped per pixel) calibration database"
0422           << std::endl;
0423     if (!theGainCalibrationDbInputHLT.put(detid, hltrange, ncols))
0424       edm::LogError("SiPixelGainCalibrationAnalysis")
0425           << "warning: detid already exists for HLT (pedestal and gain per column) calibration database" << std::endl;
0426 
0427     //delete gainPerDetid;
0428   }
0429   // now printing out summary:
0430   size_t nempty = 0;
0431   size_t ndefault = 0;
0432   size_t ndead = 0;
0433   size_t ncoldefault = 0;
0434   size_t ncoldead = 0;
0435   for (std::map<uint32_t, std::pair<TString, int> >::const_iterator ibad = badresults.begin(); ibad != badresults.end();
0436        ++ibad) {
0437     uint32_t detid = ibad->first;
0438     if (badresults[detid].second == 0) {
0439       //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << " used pixel mean value";
0440       nempty++;
0441     } else if (badresults[detid].second == 1) {
0442       //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << " used pixel mean value";
0443       ndefault++;
0444     } else if (badresults[detid].second == 2) {
0445       edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << badresults[detid].first;
0446       edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << " has one or more dead pixels";
0447       edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << std::endl;
0448       ndead++;
0449     } else if (badresults[detid].second == 3) {
0450       //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << " used column mean value";
0451       ncoldefault++;
0452     } else if (badresults[detid].second == 4) {
0453       edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << badresults[detid].first;
0454       edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << " has one or more dead columns";
0455       ncoldead++;
0456       edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << std::endl;
0457     }
0458   }
0459   edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0460       << nempty << " modules were empty and now have pixels filled with default values." << std::endl;
0461   edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0462       << ndefault << " modules have pixels filled with default values." << std::endl;
0463   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << ndead << " modules have pixels flagged as dead." << std::endl;
0464   edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0465       << ncoldefault << " modules have columns filled with default values." << std::endl;
0466   edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0467       << ncoldead << " modules have columns filled with dead values." << std::endl;
0468   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << " ---> PIXEL Modules  " << nmodules << "\n"
0469                                                      << " ---> PIXEL Channels " << nchannels << std::endl;
0470 
0471   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << " --- writing to DB!" << std::endl;
0472   edm::Service<cond::service::PoolDBOutputService> mydbservice;
0473   if (!mydbservice.isAvailable()) {
0474     edm::LogError("db service unavailable");
0475     return;
0476   } else {
0477     if (record_ == "SiPixelGainCalibrationForHLTRcd") {
0478       edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0479           << "now doing SiPixelGainCalibrationForHLTRcd payload..." << std::endl;
0480       if (mydbservice->isNewTagRequest(record_)) {
0481         mydbservice->createOneIOV<SiPixelGainCalibrationForHLT>(
0482             theGainCalibrationDbInputHLT, mydbservice->beginOfTime(), "SiPixelGainCalibrationForHLTRcd");
0483       } else {
0484         mydbservice->appendOneIOV<SiPixelGainCalibrationForHLT>(
0485             theGainCalibrationDbInputHLT, mydbservice->currentTime(), "SiPixelGainCalibrationForHLTRcd");
0486       }
0487     } else if (record_ == "SiPixelGainCalibrationOfflineRcd") {
0488       edm::LogPrint("SiPixelGainCalibrationReadDQMFile")
0489           << "now doing SiPixelGainCalibrationOfflineRcd payload..." << std::endl;
0490       if (mydbservice->isNewTagRequest(record_)) {
0491         mydbservice->createOneIOV<SiPixelGainCalibrationOffline>(
0492             theGainCalibrationDbInputOffline, mydbservice->beginOfTime(), "SiPixelGainCalibrationOfflineRcd");
0493       } else {
0494         mydbservice->appendOneIOV<SiPixelGainCalibrationOffline>(
0495             theGainCalibrationDbInputOffline, mydbservice->currentTime(), "SiPixelGainCalibrationOfflineRcd");
0496       }
0497     }
0498     edm::LogInfo(" --- all OK");
0499   }
0500 }
0501 
0502 SiPixelGainCalibrationReadDQMFile::SiPixelGainCalibrationReadDQMFile(const edm::ParameterSet &iConfig)
0503     : appendMode_(iConfig.getUntrackedParameter<bool>("appendMode", true)),
0504       pddToken_(esConsumes()),
0505       theGainCalibrationDbInputService_(iConfig, consumesCollector()),
0506       record_(iConfig.getUntrackedParameter<std::string>("record", "SiPixelGainCalibrationOfflineRcd")),
0507       gainlow_(10.),
0508       gainhi_(0.),
0509       pedlow_(255.),
0510       pedhi_(-256),
0511       usemeanwhenempty_(iConfig.getUntrackedParameter<bool>("useMeanWhenEmpty", false)),
0512       rootfilestring_(iConfig.getUntrackedParameter<std::string>("inputrootfile", "inputfile.root")),
0513       gainmax_(20),
0514       pedmax_(200),
0515       badchi2_(iConfig.getUntrackedParameter<double>("badChi2Prob", 0.01)),
0516       nmaxcols(10 * 52),
0517       nmaxrows(160) {
0518   usesResource(TFileService::kSharedResource);
0519 
0520   //now do what ever initialization is needed
0521   ::putenv((char *)"CORAL_AUTH_USER=me");
0522   ::putenv((char *)"CORAL_AUTH_PASSWORD=test");
0523   meanGainHist_ = std::make_unique<TH1F>("meanGainHist", "mean Gain Hist", 500, 0, gainmax_);
0524   meanPedHist_ = std::make_unique<TH1F>("meanPedHist", "mean Ped Hist", 512, -200, pedmax_);
0525   defaultGain_ = std::make_unique<TH2F>("defaultGain",
0526                                         "default gain, contains mean",
0527                                         nmaxcols,
0528                                         0,
0529                                         nmaxcols,
0530                                         nmaxrows,
0531                                         0,
0532                                         nmaxrows);  // using dummy (largest) module size
0533   defaultPed_ = std::make_unique<TH2F>("defaultPed",
0534                                        "default pedestal, contains mean",
0535                                        nmaxcols,
0536                                        0,
0537                                        nmaxcols,
0538                                        nmaxrows,
0539                                        0,
0540                                        nmaxrows);  // using dummy (largest) module size
0541   defaultFitResult_ = std::make_unique<TH2F>("defaultFitResult",
0542                                              "default fitresult, contains '0'",
0543                                              nmaxcols,
0544                                              0,
0545                                              nmaxcols,
0546                                              nmaxrows,
0547                                              0,
0548                                              nmaxrows);  // using dummy (largest) module size
0549   defaultChi2_ = std::make_unique<TH2F>("defaultChi2",
0550                                         "default chi2 probability, contains '1'",
0551                                         nmaxcols,
0552                                         0,
0553                                         nmaxcols,
0554                                         nmaxrows,
0555                                         0,
0556                                         nmaxrows);  // using dummy (largest) module size
0557 }
0558 
0559 //
0560 // member functions
0561 //
0562 
0563 // ------------ method called to for each event  ------------
0564 void SiPixelGainCalibrationReadDQMFile::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0565   auto histfile = getHistograms();
0566   fillDatabase(iSetup, histfile.get());
0567   // empty but should not be called anyway
0568 }
0569 
0570 std::unique_ptr<TFile> SiPixelGainCalibrationReadDQMFile::getHistograms() {
0571   edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now parsing file " << rootfilestring_ << std::endl;
0572   auto therootfile = std::make_unique<TFile>(rootfilestring_.c_str());
0573   therootfile->cd();
0574   TDirectory *dir = therootfile->GetDirectory("DQMData");
0575   TList *list = dir->GetListOfKeys();
0576 
0577   TString comparestring = "Module";
0578 
0579   std::vector<TString> keylist;
0580   std::vector<TString> hist2list;
0581   std::vector<TString> dirlist;
0582   std::vector<TString> notdonelist;
0583   std::vector<int> nsubdirs;
0584   int ikey = 0;
0585 
0586   for (ikey = 0; ikey < list->GetEntries(); ikey++) {
0587     TKey *thekey = (TKey *)list->At(ikey);
0588     if (thekey == nullptr)
0589       continue;
0590     TString keyname = thekey->GetName();
0591     TString keytype = thekey->GetClassName();
0592     //    if(keyname=="EventInfo")
0593     //      continue;
0594     //    edm::LogPrint("SiPixelGainCalibrationReadDQMFile") <<  keytype << " " << keyname << std::endl;
0595     if (keytype == "TDirectoryFile") {
0596       TString dirname = dir->GetPath();
0597       dirname += "/";
0598       dirname += keyname;
0599       //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << dirname << std::endl;
0600       dir = therootfile->GetDirectory(dirname);
0601 
0602       list = dir->GetListOfKeys();
0603       if (dirname.Contains(comparestring)) {
0604         dirlist.push_back(dirname);
0605         //  edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << dirname << std::endl;
0606       } else {
0607         notdonelist.push_back(dirname);
0608         nsubdirs.push_back(-1);
0609       }
0610     }
0611   }
0612   size_t nempty = 0;
0613   while (nempty != notdonelist.size()) {
0614     for (size_t idir = 0; idir < notdonelist.size(); ++idir) {
0615       if (nsubdirs[idir] == 0)
0616         continue;
0617       //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now examining " << notdonelist[idir]<< " " << nsubdirs[idir] <<  std::endl;
0618       dir = therootfile->GetDirectory(notdonelist[idir]);
0619       //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << dir->GetName() << std::endl;
0620       list = dir->GetListOfKeys();
0621       //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << list->GetEntries() << std::endl;
0622       int ndirectories = 0;
0623       for (ikey = 0; ikey < list->GetEntries(); ikey++) {
0624         TKey *thekey = (TKey *)list->At(ikey);
0625         if (thekey == nullptr)
0626           continue;
0627         TString keyname = thekey->GetName();
0628         TString keytype = thekey->GetClassName();
0629         //  edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << keyname << " " << keytype << std::endl;
0630         if (keytype == "TDirectoryFile") {
0631           TString dirname = dir->GetPath();
0632           dirname += "/";
0633           dirname += keyname;
0634           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << dirname << std::endl;
0635           ndirectories++;
0636           if (dirname.Contains(comparestring)) {
0637             //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << dirname << std::endl;
0638             dirlist.push_back(dirname);
0639           } else {
0640             notdonelist.push_back(dirname);
0641             nsubdirs.push_back(-1);
0642           }
0643         }
0644       }
0645       nsubdirs[idir] = ndirectories;
0646       // edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "now done examining " << notdonelist[idir]<< " " << nsubdirs[idir] <<  std::endl;
0647     }
0648     nempty = 0;
0649     for (size_t i = 0; i < nsubdirs.size(); i++) {
0650       if (nsubdirs[i] != -1)
0651         nempty++;
0652     }
0653   }
0654   //  edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "\n done!" << std::endl;
0655 
0656   for (size_t idir = 0; idir < dirlist.size(); ++idir) {
0657     //    edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << "good dir "  << dirlist[idir] << std::endl;
0658 
0659     uint32_t detid = 1;
0660 
0661     dir = therootfile->GetDirectory(dirlist[idir]);
0662     list = dir->GetListOfKeys();
0663     for (ikey = 0; ikey < list->GetEntries(); ikey++) {
0664       TKey *thekey = (TKey *)list->At(ikey);
0665       if (thekey == nullptr)
0666         continue;
0667       TString keyname = thekey->GetName();
0668       TString keytype = thekey->GetClassName();
0669       //    edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << keyname << " " << keytype << std::endl;
0670       if (keytype == "TH2F" && (keyname.Contains("Gain2d") || keyname.Contains("Pedestal2d") ||
0671                                 keyname.Contains("GainChi2Prob2d") || keyname.Contains("GainFitResult2d"))) {
0672         TString detidstring = keyname;
0673         detidstring.Remove(0, detidstring.Sizeof() - 10);
0674 
0675         detid = atoi(detidstring.Data());
0676 
0677         if (keyname.Contains("GainChi2Prob2d")) {
0678           TString tempstr = dirlist[idir];
0679           tempstr += "/";
0680           tempstr += keyname;
0681           TString replacestring = rootfilestring_;
0682           replacestring += ":";
0683           tempstr.ReplaceAll(replacestring, "");
0684           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << tempstr << std::endl;
0685           bookkeeper_[detid]["chi2prob_2d"] = tempstr;
0686         } else if (keyname.Contains("GainFitResult2d")) {
0687           TString tempstr = dirlist[idir];
0688           tempstr += "/";
0689           tempstr += keyname;
0690           TString replacestring = rootfilestring_;
0691           replacestring += ":";
0692           tempstr.ReplaceAll(replacestring, "");
0693           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << tempstr << std::endl;
0694           bookkeeper_[detid]["fitresult_2d"] = tempstr;
0695         } else if (keyname.Contains("Gain2d")) {
0696           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << dirlist[idir] << std::endl;
0697           std::map<std::string, TString> tempmap;
0698           TString tempstr = dirlist[idir];
0699           tempstr += "/";
0700           tempstr += keyname;
0701           TString replacestring = rootfilestring_;
0702           replacestring += ":";
0703           tempstr.ReplaceAll(replacestring, "");
0704           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << tempstr << std::endl;
0705           bookkeeper_[detid]["gain_2d"] = tempstr;
0706           //edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << detidstring << " " << keyname << " " << detid << " " << bookkeeper_[detid]["gain_2d"] << std::endl ;
0707         }
0708         if (keyname.Contains("Pedestal2d")) {
0709           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << dirlist[idir] << std::endl;
0710           std::map<std::string, TString> tempmap;
0711 
0712           TString tempstr = dirlist[idir];
0713           tempstr += "/";
0714           tempstr += keyname;
0715           //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << tempstr << std::endl;
0716           TString replacestring = rootfilestring_;
0717           replacestring += ":";
0718           tempstr.ReplaceAll(replacestring, "");
0719           bookkeeper_[detid]["ped_2d"] = tempstr;
0720 
0721           //edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << detidstring << " " << keyname << " " << detid << " " << bookkeeper_[detid]["ped_2d"] << std::endl ;
0722         }
0723         //      edm::LogPrint("SiPixelGainCalibrationReadDQMFile") << keyname << " " << keytype << std::endl;
0724       }
0725     }
0726 
0727     TH2F *temphistoped = dynamic_cast<TH2F *>(therootfile->Get(bookkeeper_[detid]["ped_2d"]));
0728     TH2F *temphistogain = dynamic_cast<TH2F *>(therootfile->Get(bookkeeper_[detid]["gain_2d"]));
0729     TH2F *temphistofitresult = dynamic_cast<TH2F *>(therootfile->Get(bookkeeper_[detid]["gain_2d"]));
0730 
0731     for (int xbin = 1; xbin <= temphistoped->GetNbinsX(); ++xbin) {
0732       for (int ybin = 1; ybin <= temphistoped->GetNbinsY(); ++ybin) {
0733         if (temphistofitresult->GetBinContent(xbin, ybin) <= 0)
0734           continue;
0735         float val = temphistoped->GetBinContent(xbin, ybin);
0736         if (val > pedmax_)
0737           continue;
0738         if (pedlow_ > val)
0739           pedlow_ = val;
0740         if (pedhi_ < val)
0741           pedhi_ = val;
0742         meanPedHist_->Fill(val);
0743       }
0744     }
0745 
0746     for (int xbin = 1; xbin <= temphistogain->GetNbinsX(); ++xbin) {
0747       for (int ybin = 1; ybin <= temphistogain->GetNbinsY(); ++ybin) {
0748         if (temphistofitresult->GetBinContent(xbin, ybin) <= 0)
0749           continue;
0750         float val = temphistogain->GetBinContent(xbin, ybin);
0751         if (val <= 0.0001)
0752           continue;
0753         if (gainlow_ > val)
0754           gainlow_ = val;
0755         if (gainhi_ < val)
0756           gainhi_ = val;
0757         meanGainHist_->Fill(val);
0758       }
0759     }
0760 
0761   }  // end of loop over dirlist
0762   return therootfile;
0763 }
0764 
0765 //define this as a plug-in
0766 DEFINE_FWK_MODULE(SiPixelGainCalibrationReadDQMFile);