Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:48:45

0001 // system includes
0002 #include <fstream>
0003 #include <iostream>
0004 #include <map>
0005 #include <memory>
0006 #include <string>
0007 
0008 // user includes
0009 #include "CondFormats/SiPixelObjects/interface/SiPixelVCal.h"
0010 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0011 #include "DataFormats/DetId/interface/DetId.h"
0012 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0013 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0014 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/Framework/interface/MakerMacros.h"
0018 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/ServiceRegistry/interface/Service.h"
0022 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0023 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0024 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0025 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0026 
0027 class SiPixelVCalDB : public edm::one::EDAnalyzer<> {
0028 public:
0029   explicit SiPixelVCalDB(const edm::ParameterSet& conf);
0030   explicit SiPixelVCalDB();
0031   ~SiPixelVCalDB() override;
0032   void analyze(const edm::Event&, const edm::EventSetup&) override;
0033 
0034 private:
0035   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
0036   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tkTopoToken_;
0037   std::string recordName_;
0038   typedef std::vector<edm::ParameterSet> Parameters;
0039   Parameters BPixParameters_;
0040   Parameters FPixParameters_;
0041 };
0042 
0043 using namespace std;
0044 using namespace edm;
0045 
0046 SiPixelVCalDB::SiPixelVCalDB(edm::ParameterSet const& iConfig)
0047     : tkGeomToken_(esConsumes()), tkTopoToken_(esConsumes()) {
0048   recordName_ = iConfig.getUntrackedParameter<std::string>("record", "SiPixelVCalRcd");
0049   BPixParameters_ = iConfig.getUntrackedParameter<Parameters>("BPixParameters");
0050   FPixParameters_ = iConfig.getUntrackedParameter<Parameters>("FPixParameters");
0051 }
0052 
0053 SiPixelVCalDB::~SiPixelVCalDB() = default;
0054 
0055 // Analyzer: Functions that gets called by framework every event
0056 void SiPixelVCalDB::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
0057   SiPixelVCal vcal;
0058   bool phase1 = true;
0059 
0060   // Retrieve tracker topology from geometry
0061   const TrackerTopology* const tTopo = &iSetup.getData(tkTopoToken_);
0062 
0063   // Retrieve old style tracker geometry from geometry
0064   const TrackerGeometry* pDD = &iSetup.getData(tkGeomToken_);
0065   edm::LogPrint("SiPixelVCalDB") << " There are " << pDD->detUnits().size() << " modules" << std::endl;
0066 
0067   for (const auto& it : pDD->detUnits()) {
0068     if (dynamic_cast<PixelGeomDetUnit const*>(it) != nullptr) {
0069       const DetId detid = it->geographicalId();
0070       const unsigned int rawDetId = detid.rawId();
0071       int subid = detid.subdetId();
0072 
0073       // FILL BPIX
0074       if (subid == static_cast<int>(PixelSubdetector::PixelBarrel)) {
0075         int layer = tTopo->pxbLayer(detid);    // 1, 2, 3, 4
0076         int ladder = tTopo->pxbLadder(detid);  // 1-12/28/44/64
0077         edm::LogPrint("SiPixelVCalDB") << " pixel barrel:"
0078                                        << " detId=" << rawDetId << ", layer=" << layer << ", ladder=" << ladder;
0079         for (Parameters::iterator it = BPixParameters_.begin(); it != BPixParameters_.end(); ++it) {
0080           if (it->getParameter<int>("layer") == layer && it->getParameter<int>("ladder") == ladder) {
0081             float slope = (float)it->getParameter<double>("slope");
0082             float offset = (float)it->getParameter<double>("offset");
0083             edm::LogPrint("SiPixelVCalDB") << ";  VCal slope " << slope << ", offset " << offset;
0084             // edm::LogInfo("SiPixelVCalDB")  << "  detId " << rawDetId << " \t
0085             // VCal slope " << slope << ", offset " << offset;
0086             vcal.putSlopeAndOffset(detid, slope, offset);
0087           }
0088         }
0089         edm::LogPrint("SiPixelVCalDB") << std::endl;
0090 
0091         // FILL FPIX
0092       } else if (subid == static_cast<int>(PixelSubdetector::PixelEndcap)) {
0093         PixelEndcapName fpix(detid, tTopo, phase1);
0094         int side = tTopo->pxfSide(detid);   // 1 (-z), 2 for (+z)
0095         int disk = fpix.diskName();         // 1, 2, 3
0096         int disk2 = tTopo->pxfDisk(detid);  // 1, 2, 3
0097         int ring = fpix.ringName();         // 1 (lower), 2 (upper)
0098         if (disk != disk2) {
0099           edm::LogError("SiPixelVCalDB::analyze")
0100               << "Found contradicting FPIX disk number: " << disk << " vs." << disk2 << std::endl;
0101         }
0102         edm::LogPrint("SiPixelVCalDB") << " pixel endcap:"
0103                                        << " detId=" << rawDetId << ", side=" << side << ", disk=" << disk
0104                                        << ", ring=" << ring;
0105         for (Parameters::iterator it = FPixParameters_.begin(); it != FPixParameters_.end(); ++it) {
0106           if (it->getParameter<int>("side") == side && it->getParameter<int>("disk") == disk &&
0107               it->getParameter<int>("ring") == ring) {
0108             float slope = (float)it->getParameter<double>("slope");
0109             float offset = (float)it->getParameter<double>("offset");
0110             edm::LogPrint("SiPixelVCalDB") << ";  VCal slope " << slope << ", offset " << offset;
0111             // edm::LogInfo("SiPixelVCalDB")  << "  detId " << rawDetId << " \t
0112             // VCal slope " << slope << ", offset " << offset;
0113             vcal.putSlopeAndOffset(rawDetId, slope, offset);
0114           }
0115         }
0116         edm::LogPrint("SiPixelVCalDB") << std::endl;
0117 
0118       } else {
0119         edm::LogError("SiPixelVCalDB::analyze") << "detid is Pixel but neither bpix nor fpix" << std::endl;
0120       }
0121     }
0122   }
0123 
0124   // Save to DB
0125   edm::Service<cond::service::PoolDBOutputService> mydbservice;
0126   if (mydbservice.isAvailable()) {
0127     try {
0128       if (mydbservice->isNewTagRequest(recordName_)) {
0129         mydbservice->createOneIOV<SiPixelVCal>(vcal, mydbservice->beginOfTime(), recordName_);
0130       } else {
0131         mydbservice->appendOneIOV<SiPixelVCal>(vcal, mydbservice->currentTime(), recordName_);
0132       }
0133     } catch (const cond::Exception& er) {
0134       edm::LogError("SiPixelVCalDB") << er.what() << std::endl;
0135     } catch (const std::exception& er) {
0136       edm::LogError("SiPixelVCalDB") << "caught std::exception " << er.what() << std::endl;
0137     } catch (...) {
0138       edm::LogError("SiPixelVCalDB") << "Funny error" << std::endl;
0139     }
0140   } else {
0141     edm::LogError("SiPixelVCalDB") << "Service is unavailable" << std::endl;
0142   }
0143 }
0144 DEFINE_FWK_MODULE(SiPixelVCalDB);