Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-17 04:57:59

0001 // -*- C++ -*-
0002 //
0003 // Package:    CalibTracker/SiPixelGainForHLTESProducers
0004 // Class:      SiPixelFakeGainForHLTESSource
0005 //
0006 /**\class SiPixelFakeGainForHLTESSource SiPixelFakeGainForHLTESSource.cc CalibTracker/SiPixelGainForHLTESProducers/plugins/SiPixelFakeGainForHLTESSource.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Vincenzo Chiochia
0015 //         Created:  Tue 8 12:31:25 CEST 2007
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"
0024 #include "CondFormats/DataRecord/interface/SiPixelGainCalibrationForHLTRcd.h"
0025 #include "CondFormats/SiPixelObjects/interface/SiPixelGainCalibrationForHLT.h"
0026 #include "FWCore/Framework/interface/ESProducer.h"
0027 #include "FWCore/Framework/interface/EventSetupRecordIntervalFinder.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/Framework/interface/ModuleFactory.h"
0030 #include "FWCore/Framework/interface/SourceFactory.h"
0031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0032 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0035 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0036 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0037 
0038 //
0039 // class decleration
0040 //
0041 
0042 class SiPixelFakeGainForHLTESSource : public edm::ESProducer, public edm::EventSetupRecordIntervalFinder {
0043 public:
0044   SiPixelFakeGainForHLTESSource(const edm::ParameterSet&);
0045   ~SiPixelFakeGainForHLTESSource() override = default;
0046 
0047   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0048 
0049   virtual std::unique_ptr<SiPixelGainCalibrationForHLT> produce(const SiPixelGainCalibrationForHLTRcd&);
0050 
0051 protected:
0052   void setIntervalFor(const edm::eventsetup::EventSetupRecordKey&,
0053                       const edm::IOVSyncValue&,
0054                       edm::ValidityInterval&) override;
0055 
0056 private:
0057   edm::FileInPath fp_;
0058 };
0059 
0060 //
0061 // constructors and destructor
0062 //
0063 SiPixelFakeGainForHLTESSource::SiPixelFakeGainForHLTESSource(const edm::ParameterSet& conf_)
0064     : fp_(conf_.getParameter<edm::FileInPath>("file")) {
0065   edm::LogInfo("SiPixelFakeGainForHLTESSource::SiPixelFakeGainForHLTESSource");
0066   //the following line is needed to tell the framework what
0067   // data is being produced
0068   setWhatProduced(this);
0069   findingRecord<SiPixelGainCalibrationForHLTRcd>();
0070 }
0071 
0072 std::unique_ptr<SiPixelGainCalibrationForHLT> SiPixelFakeGainForHLTESSource::produce(
0073     const SiPixelGainCalibrationForHLTRcd&) {
0074   using namespace edm::es;
0075   SiPixelGainCalibrationForHLT* obj = new SiPixelGainCalibrationForHLT(25., 30., 2., 3.);
0076   SiPixelDetInfoFileReader reader(fp_.fullPath());
0077   const std::vector<uint32_t>& DetIds = reader.getAllDetIds();
0078 
0079   // Loop over detectors
0080   for (std::vector<uint32_t>::const_iterator detit = DetIds.begin(); detit != DetIds.end(); detit++) {
0081     std::vector<char> theSiPixelGainCalibration;
0082     const std::pair<int, int>& detUnitDimensions = reader.getDetUnitDimensions(*detit);
0083 
0084     // Loop over columns and rows
0085 
0086     for (int i = 0; i < detUnitDimensions.first; i++) {
0087       float totalGain = 0.0;
0088       float totalPed = 0.0;
0089       float totalEntries = 0.0;
0090       for (int j = 0; j < detUnitDimensions.second; j++) {
0091         totalGain += 2.8;
0092         totalPed += 28.2;
0093         totalEntries++;
0094         if ((j + 1) % 80 == 0) {
0095           float gain = totalGain / totalEntries;
0096           float ped = totalPed / totalEntries;
0097 
0098           obj->setData(ped, gain, theSiPixelGainCalibration);
0099           totalGain = 0.;
0100           totalPed = 0.;
0101           totalEntries = 0.;
0102         }
0103       }
0104     }
0105 
0106     //std::cout << "detid " << (*detit) << std::endl;
0107 
0108     SiPixelGainCalibrationForHLT::Range range(theSiPixelGainCalibration.begin(), theSiPixelGainCalibration.end());
0109     int nCols = detUnitDimensions.first;
0110     if (!obj->put(*detit, range, nCols))
0111       edm::LogError("SiPixelFakeGainForHLTESSource")
0112           << "[SiPixelFakeGainForHLTESSource::produce] detid already exists" << std::endl;
0113   }
0114 
0115   //
0116   return std::unique_ptr<SiPixelGainCalibrationForHLT>(obj);
0117 }
0118 
0119 void SiPixelFakeGainForHLTESSource::setIntervalFor(const edm::eventsetup::EventSetupRecordKey&,
0120                                                    const edm::IOVSyncValue& iosv,
0121                                                    edm::ValidityInterval& oValidity) {
0122   edm::ValidityInterval infinity(iosv.beginOfTime(), iosv.endOfTime());
0123   oValidity = infinity;
0124 }
0125 
0126 void SiPixelFakeGainForHLTESSource::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0127   edm::ParameterSetDescription desc;
0128   desc.add<edm::FileInPath>("file", edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt"));
0129   descriptions.addWithDefaultLabel(desc);
0130 }
0131 
0132 DEFINE_FWK_EVENTSETUP_SOURCE(SiPixelFakeGainForHLTESSource);