Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiPixelCondObjForHLTBuilder
0004 // Class:      SiPixelCondObjForHLTBuilder
0005 //
0006 /**\class SiPixelCondObjForHLTBuilder SiPixelCondObjForHLTBuilder.h SiPixel/test/SiPixelCondObjForHLTBuilder.h
0007 
0008  Description: Test analyzer for writing pixel calibration in the DB
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Vincenzo CHIOCHIA
0015 //         Created:  Tue Oct 17 17:40:56 CEST 2006
0016 // $Id: SiPixelCondObjForHLTBuilder.h,v 1.6 2009/11/20 19:21:02 rougny Exp $
0017 //
0018 //
0019 
0020 // system includes
0021 #include <string>
0022 #include <memory>
0023 #include <iostream>
0024 
0025 // user includes
0026 #include "CalibTracker/SiPixelESProducers/interface/SiPixelGainCalibrationForHLTService.h"
0027 #include "CondFormats/SiPixelObjects/interface/PixelIndices.h"
0028 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0029 #include "DataFormats/DetId/interface/DetId.h"
0030 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0031 #include "FWCore/Framework/interface/Event.h"
0032 #include "FWCore/Framework/interface/EventSetup.h"
0033 #include "FWCore/Framework/interface/Frameworkfwd.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0038 #include "FWCore/ServiceRegistry/interface/Service.h"
0039 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0040 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0041 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0042 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0043 
0044 #include "CLHEP/Random/RandGauss.h"
0045 #include "CLHEP/Random/RandFlat.h"
0046 
0047 namespace cms {
0048   class SiPixelCondObjForHLTBuilder : public edm::one::EDAnalyzer<> {
0049   public:
0050     explicit SiPixelCondObjForHLTBuilder(const edm::ParameterSet& iConfig);
0051     ~SiPixelCondObjForHLTBuilder() override = default;
0052     void beginJob() override;
0053     void analyze(const edm::Event&, const edm::EventSetup&) override;
0054     bool loadFromFile();
0055 
0056   private:
0057     const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeometryToken_;
0058 
0059     bool appendMode_;
0060     std::unique_ptr<SiPixelGainCalibrationForHLT> SiPixelGainCalibration_;
0061     SiPixelGainCalibrationForHLTService SiPixelGainCalibrationService_;
0062     const std::string recordName_;
0063 
0064     const double meanPed_;
0065     const double rmsPed_;
0066     const double meanGain_;
0067     const double rmsGain_;
0068     const double meanPedFPix_;
0069     const double rmsPedFPix_;
0070     const double meanGainFPix_;
0071     const double rmsGainFPix_;
0072     const double deadFraction_;
0073     const double noisyFraction_;
0074     const double secondRocRowGainOffset_;
0075     const double secondRocRowPedOffset_;
0076     const int numberOfModules_;
0077     const bool fromFile_;
0078     const std::string fileName_;
0079     const bool generateColumns_;
0080 
0081     // Internal class
0082     class CalParameters {
0083     public:
0084       float p0;
0085       float p1;
0086     };
0087     // Map for storing calibration constants
0088     std::map<int, CalParameters, std::less<int> > calmap_;
0089     PixelIndices* pIndexConverter_;  // Pointer to the index converter
0090   };
0091 }  // namespace cms
0092 
0093 namespace cms {
0094   SiPixelCondObjForHLTBuilder::SiPixelCondObjForHLTBuilder(const edm::ParameterSet& iConfig)
0095       : tkGeometryToken_(esConsumes()),
0096         appendMode_(iConfig.getUntrackedParameter<bool>("appendMode", true)),
0097         SiPixelGainCalibration_(nullptr),
0098         SiPixelGainCalibrationService_(iConfig, consumesCollector()),
0099         recordName_(iConfig.getParameter<std::string>("record")),
0100         meanPed_(iConfig.getParameter<double>("meanPed")),
0101         rmsPed_(iConfig.getParameter<double>("rmsPed")),
0102         meanGain_(iConfig.getParameter<double>("meanGain")),
0103         rmsGain_(iConfig.getParameter<double>("rmsGain")),
0104         meanPedFPix_(iConfig.getUntrackedParameter<double>("meanPedFPix", meanPed_)),
0105         rmsPedFPix_(iConfig.getUntrackedParameter<double>("rmsPedFPix", rmsPed_)),
0106         meanGainFPix_(iConfig.getUntrackedParameter<double>("meanGainFPix", meanGain_)),
0107         rmsGainFPix_(iConfig.getUntrackedParameter<double>("rmsGainFPix", rmsGain_)),
0108         deadFraction_(iConfig.getParameter<double>("deadFraction")),
0109         noisyFraction_(iConfig.getParameter<double>("noisyFraction")),
0110         secondRocRowGainOffset_(iConfig.getParameter<double>("secondRocRowGainOffset")),
0111         secondRocRowPedOffset_(iConfig.getParameter<double>("secondRocRowPedOffset")),
0112         numberOfModules_(iConfig.getParameter<int>("numberOfModules")),
0113         fromFile_(iConfig.getParameter<bool>("fromFile")),
0114         fileName_(iConfig.getParameter<std::string>("fileName")),
0115         generateColumns_(iConfig.getUntrackedParameter<bool>("generateColumns", true))
0116 
0117   {
0118     ::putenv((char*)"CORAL_AUTH_USER=me");
0119     ::putenv((char*)"CORAL_AUTH_PASSWORD=test");
0120   }
0121 
0122   void SiPixelCondObjForHLTBuilder::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0123     using namespace edm;
0124     unsigned int run = iEvent.id().run();
0125     int nmodules = 0;
0126     uint32_t nchannels = 0;
0127     //    int mycol = 415;
0128     //    int myrow = 159;
0129 
0130     edm::LogInfo("SiPixelCondObjForHLTBuilder")
0131         << "... creating dummy SiPixelGainCalibration Data for Run " << run << "\n " << std::endl;
0132     //
0133     // Instantiate Gain calibration offset and define pedestal/gain range
0134     //
0135     // note: the hard-coded range values are also used in random generation. That is why they're defined here
0136 
0137     float mingain = 0;
0138     float maxgain = 10;
0139     float minped = 0;
0140     float maxped = 255;
0141     SiPixelGainCalibration_ = std::make_unique<SiPixelGainCalibrationForHLT>(minped, maxped, mingain, maxgain);
0142 
0143     const TrackerGeometry* pDD = &iSetup.getData(tkGeometryToken_);
0144     edm::LogInfo("SiPixelCondObjForHLTBuilder") << " There are " << pDD->dets().size() << " detectors" << std::endl;
0145 
0146     for (TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++) {
0147       if (dynamic_cast<PixelGeomDetUnit const*>((*it)) != nullptr) {
0148         uint32_t detid = ((*it)->geographicalId()).rawId();
0149 
0150         // Stop if module limit reached
0151         nmodules++;
0152         if (nmodules > numberOfModules_)
0153           break;
0154 
0155         const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>((*it));
0156         const PixelTopology& topol = pixDet->specificTopology();
0157         // Get the module sizes.
0158         int nrows = topol.nrows();     // rows in x
0159         int ncols = topol.ncolumns();  // cols in y
0160         //edm::LogPrint("SiPixelCondObjForHLTBuilder") << " ---> PIXEL DETID " << detid << " Cols " << ncols << " Rows " << nrows << std::endl;
0161 
0162         double meanPedWork = meanPed_;
0163         double rmsPedWork = rmsPed_;
0164         double meanGainWork = meanGain_;
0165         double rmsGainWork = rmsGain_;
0166         DetId detId(detid);
0167         if (detId.subdetId() == 2) {  // FPIX
0168           meanPedWork = meanPedFPix_;
0169           rmsPedWork = rmsPedFPix_;
0170           meanGainWork = meanGainFPix_;
0171           rmsGainWork = rmsGainFPix_;
0172         }
0173 
0174         PixelIndices pIndexConverter(ncols, nrows);
0175 
0176         std::vector<char> theSiPixelGainCalibration;
0177 
0178         // Loop over columns and rows of this DetID
0179         for (int i = 0; i < ncols; i++) {
0180           float totalPed = 0.0;
0181           float totalGain = 0.0;
0182           for (int j = 0; j < nrows; j++) {
0183             nchannels++;
0184             bool isDead = false;
0185             bool isNoisy = false;
0186             float ped = 0.0, gain = 0.0;
0187 
0188             if (fromFile_) {
0189               // Use calibration from a file
0190               int chipIndex = 0, colROC = 0, rowROC = 0;
0191 
0192               pIndexConverter.transformToROC(i, j, chipIndex, colROC, rowROC);
0193               int chanROC = PixelIndices::pixelToChannelROC(rowROC, colROC);  // use ROC coordinates
0194               //         float pp0=0, pp1=0;
0195               std::map<int, CalParameters, std::less<int> >::const_iterator it = calmap_.find(chanROC);
0196               CalParameters theCalParameters = (*it).second;
0197               ped = theCalParameters.p0;
0198               gain = theCalParameters.p1;
0199 
0200             } else {
0201               if (deadFraction_ > 0) {
0202                 double val = CLHEP::RandFlat::shoot();
0203                 if (val < deadFraction_) {
0204                   isDead = true;
0205                   //         edm::LogPrint("SiPixelCondObjForHLTBuilder") << "dead pixel " << detid << " " << i << "," << j << " " << val << std::endl;
0206                 }
0207               }
0208               if (deadFraction_ > 0 && !isDead) {
0209                 double val = CLHEP::RandFlat::shoot();
0210                 if (val < noisyFraction_) {
0211                   isNoisy = true;
0212                   //         edm::LogPrint("SiPixelCondObjForHLTBuilder") << "noisy pixel " << detid << " " << i << "," << j << " " << val << std::endl;
0213                 }
0214               }
0215 
0216               if (rmsPedWork > 0) {
0217                 ped = CLHEP::RandGauss::shoot(meanPedWork, rmsPedWork);
0218                 while (ped < minped || ped > maxped)
0219                   ped = CLHEP::RandGauss::shoot(meanPedWork, rmsPedWork);
0220               } else
0221                 ped = meanPedWork;
0222               if (rmsGainWork > 0) {
0223                 gain = CLHEP::RandGauss::shoot(meanGainWork, rmsGainWork);
0224                 while (gain < mingain || gain > maxgain)
0225                   gain = CLHEP::RandGauss::shoot(meanGainWork, rmsGainWork);
0226               } else
0227                 gain = meanGainWork;
0228             }
0229 
0230             //     if(i==mycol && j==myrow) {
0231             //     }
0232 
0233             //     gain =  2.8;
0234             //     ped  = 28.2;
0235 
0236             //if in the second row of rocs (i.e. a 2xN plaquette) add an offset (if desired) for testing
0237             if (j >= 80) {
0238               ped += secondRocRowPedOffset_;
0239               gain += secondRocRowGainOffset_;
0240 
0241               if (gain > maxgain)
0242                 gain = maxgain;
0243               else if (gain < mingain)
0244                 gain = mingain;
0245 
0246               if (ped > maxped)
0247                 ped = maxped;
0248               else if (ped < minped)
0249                 ped = minped;
0250             }
0251 
0252             totalPed += ped;
0253             totalGain += gain;
0254 
0255             if ((j + 1) % 80 == 0) {
0256               //edm::LogPrint("SiPixelCondObjForHLTBuilder") << "Filling   Col "<<i<<" Row "<<j<<" Ped "<<totalPed<<" Gain "<<totalGain<<std::endl;
0257               float averagePed = totalPed / static_cast<float>(80);
0258               float averageGain = totalGain / static_cast<float>(80);
0259 
0260               if (generateColumns_) {
0261                 averagePed = ped;
0262                 averageGain = gain;
0263               }
0264               //only fill by column after each roc
0265               if (!isDead && !isNoisy)
0266                 SiPixelGainCalibration_->setData(averagePed, averageGain, theSiPixelGainCalibration);
0267               else if (isDead)
0268                 SiPixelGainCalibration_->setDeadColumn(80, theSiPixelGainCalibration);
0269               else if (isNoisy)
0270                 SiPixelGainCalibration_->setNoisyColumn(80, theSiPixelGainCalibration);
0271 
0272               totalPed = 0;
0273               totalGain = 0;
0274             }
0275           }
0276         }
0277 
0278         SiPixelGainCalibrationForHLT::Range range(theSiPixelGainCalibration.begin(), theSiPixelGainCalibration.end());
0279         if (!SiPixelGainCalibration_->put(detid, range, ncols))
0280           edm::LogError("SiPixelCondObjForHLTBuilder")
0281               << "[SiPixelCondObjForHLTBuilder::analyze] detid already exists" << std::endl;
0282       }
0283     }
0284     edm::LogPrint("SiPixelCondObjForHLTBuilder") << " ---> PIXEL Modules  " << nmodules << std::endl;
0285     edm::LogPrint("SiPixelCondObjForHLTBuilder") << " ---> PIXEL Channels " << nchannels << std::endl;
0286 
0287     //   // Try to read object
0288     //    int mynmodules =0;
0289     //    for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
0290     //      if( dynamic_cast<PixelGeomDetUnit const*>((*it))!=0){
0291     //        uint32_t mydetid=((*it)->geographicalId()).rawId();
0292     //        mynmodules++;
0293     //        if( mynmodules > numberOfModules_) break;
0294     //        SiPixelGainCalibration::Range myrange = SiPixelGainCalibration_->getRange(mydetid);
0295     //        float mypedestal = SiPixelGainCalibration_->getPed (mycol,myrow,myrange,416);
0296     //        float mygain     = SiPixelGainCalibration_->getGain(mycol,myrow,myrange,416);
0297     //        //edm::LogPrint("SiPixelCondObjForHLTBuilder")<<" PEDESTAL "<< mypedestal<<" GAIN "<<mygain<<std::endl;
0298     //      }
0299     //    }
0300     // Write into DB
0301     edm::LogInfo(" --- writeing to DB!");
0302     edm::Service<cond::service::PoolDBOutputService> mydbservice;
0303     if (!mydbservice.isAvailable()) {
0304       edm::LogError("db service unavailable");
0305       return;
0306     } else {
0307       edm::LogInfo("DB service OK");
0308     }
0309 
0310     try {
0311       //     size_t callbackToken = mydbservice->callbackToken("SiPixelGainCalibration");
0312       //     edm::LogInfo("SiPixelCondObjForHLTBuilder")<<"CallbackToken SiPixelGainCalibration "
0313       //         <<callbackToken<<std::endl;
0314       //       unsigned long long tillTime;
0315       //     if ( appendMode_)
0316       //     tillTime = mydbservice->currentTime();
0317       //       else
0318       //     tillTime = mydbservice->endOfTime();
0319       //     edm::LogInfo("SiPixelCondObjForHLTBuilder")<<"[SiPixelCondObjForHLTBuilder::analyze] tillTime = "
0320       //         <<tillTime<<std::endl;
0321       //     mydbservice->newValidityForNewPayload<SiPixelGainCalibration>(
0322       //           SiPixelGainCalibration_, tillTime , callbackToken);
0323 
0324       if (mydbservice->isNewTagRequest(recordName_)) {
0325         mydbservice->createOneIOV<SiPixelGainCalibrationForHLT>(
0326             *SiPixelGainCalibration_, mydbservice->beginOfTime(), recordName_);
0327       } else {
0328         mydbservice->appendOneIOV<SiPixelGainCalibrationForHLT>(
0329             *SiPixelGainCalibration_, mydbservice->currentTime(), recordName_);
0330       }
0331       edm::LogInfo(" --- all OK");
0332     } catch (const cond::Exception& er) {
0333       edm::LogError("SiPixelCondObjForHLTBuilder") << er.what() << std::endl;
0334     } catch (const std::exception& er) {
0335       edm::LogError("SiPixelCondObjForHLTBuilder") << "caught std::exception " << er.what() << std::endl;
0336     } catch (...) {
0337       edm::LogError("SiPixelCondObjForHLTBuilder") << "Funny error" << std::endl;
0338     }
0339   }
0340 
0341   // ------------ method called once each job just before starting event loop  ------------
0342   void SiPixelCondObjForHLTBuilder::beginJob() {
0343     if (fromFile_) {
0344       if (loadFromFile()) {
0345         edm::LogInfo("SiPixelCondObjForHLTBuilder") << " Calibration loaded: Map size " << calmap_.size() << " max "
0346                                                     << calmap_.max_size() << " " << calmap_.empty() << std::endl;
0347       }
0348     }
0349   }
0350 
0351   bool SiPixelCondObjForHLTBuilder::loadFromFile() {
0352     float par0, par1;  //,par2,par3;
0353     int colid, rowid;  //rocid
0354     std::string name;
0355 
0356     std::ifstream in_file;                          // data file pointer
0357     in_file.open(fileName_.c_str(), std::ios::in);  // in C++
0358     if (in_file.bad()) {
0359       edm::LogError("SiPixelCondObjForHLTBuilder") << "Input file not found" << std::endl;
0360     }
0361     if (in_file.eof() != 0) {
0362       edm::LogError("SiPixelCondObjForHLTBuilder") << in_file.eof() << " " << in_file.gcount() << " " << in_file.fail()
0363                                                    << " " << in_file.good() << " end of file " << std::endl;
0364       return false;
0365     }
0366     //Load file header
0367     char line[500];
0368     for (int i = 0; i < 3; i++) {
0369       in_file.getline(line, 500, '\n');
0370       edm::LogInfo("SiPixelCondObjForHLTBuilder") << line << std::endl;
0371     }
0372     //Loading calibration constants from file, loop on pixels
0373     for (int i = 0; i < (52 * 80); i++) {
0374       in_file >> par0 >> par1 >> name >> colid >> rowid;
0375 
0376       edm::LogPrint("SiPixelCondObjForHLTBuilder")
0377           << " Col " << colid << " Row " << rowid << " P0 " << par0 << " P1 " << par1 << std::endl;
0378 
0379       CalParameters onePix;
0380       onePix.p0 = par0;
0381       onePix.p1 = par1;
0382 
0383       // Convert ROC pixel index to channel
0384       int chan = PixelIndices::pixelToChannelROC(rowid, colid);
0385       calmap_.insert(std::pair<int, CalParameters>(chan, onePix));
0386     }
0387 
0388     bool flag;
0389     if (calmap_.size() == 4160) {
0390       flag = true;
0391     } else {
0392       flag = false;
0393     }
0394     return flag;
0395   }
0396 
0397 }  // namespace cms
0398 
0399 using namespace cms;
0400 DEFINE_FWK_MODULE(SiPixelCondObjForHLTBuilder);