Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // system includes
0002 #include <map>
0003 #include <memory>
0004 #include <string>
0005 #include <iostream>
0006 #include <fstream>
0007 
0008 // user includes
0009 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0010 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
0011 #include "DataFormats/DetId/interface/DetId.h"
0012 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0013 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0014 #include "FWCore/Framework/interface/Event.h"
0015 #include "FWCore/Framework/interface/EventSetup.h"
0016 #include "FWCore/Framework/interface/MakerMacros.h"
0017 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 #include "FWCore/ServiceRegistry/interface/Service.h"
0021 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0022 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0023 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0024 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0025 #include "MagneticField/Engine/interface/MagneticField.h"
0026 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0027 
0028 class SiPixelLorentzAngleDB : public edm::one::EDAnalyzer<> {
0029 public:
0030   explicit SiPixelLorentzAngleDB(const edm::ParameterSet& conf);
0031   ~SiPixelLorentzAngleDB() override;
0032   void analyze(const edm::Event& e, const edm::EventSetup& c) override;
0033 
0034 private:
0035   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
0036   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tkTopoToken_;
0037 
0038   unsigned int HVgroup(unsigned int panel, unsigned int module);
0039 
0040   std::vector<std::pair<uint32_t, float> > detid_la;
0041   double magneticField_;
0042   std::string recordName_;
0043 
0044   typedef std::vector<edm::ParameterSet> Parameters;
0045   Parameters BPixParameters_;
0046   Parameters FPixParameters_;
0047   Parameters ModuleParameters_;
0048 
0049   std::string fileName_;
0050   bool useFile_;
0051 };
0052 
0053 using namespace std;
0054 using namespace edm;
0055 
0056 //Constructor
0057 
0058 SiPixelLorentzAngleDB::SiPixelLorentzAngleDB(edm::ParameterSet const& conf)
0059     : tkGeomToken_(esConsumes()), tkTopoToken_(esConsumes()) {
0060   magneticField_ = conf.getParameter<double>("magneticField");
0061   recordName_ = conf.getUntrackedParameter<std::string>("record", "SiPixelLorentzAngleRcd");
0062   useFile_ = conf.getParameter<bool>("useFile");
0063   fileName_ = conf.getParameter<string>("fileName");
0064 
0065   BPixParameters_ = conf.getUntrackedParameter<Parameters>("BPixParameters");
0066   FPixParameters_ = conf.getUntrackedParameter<Parameters>("FPixParameters");
0067   ModuleParameters_ = conf.getUntrackedParameter<Parameters>("ModuleParameters");
0068 }
0069 
0070 // Virtual destructor needed.
0071 SiPixelLorentzAngleDB::~SiPixelLorentzAngleDB() = default;
0072 
0073 // Analyzer: Functions that gets called by framework every event
0074 
0075 void SiPixelLorentzAngleDB::analyze(const edm::Event& e, const edm::EventSetup& es) {
0076   SiPixelLorentzAngle LorentzAngle;
0077 
0078   //Retrieve tracker topology from geometry
0079   const TrackerTopology* tTopo = &es.getData(tkTopoToken_);
0080 
0081   //Retrieve old style tracker geometry from geometry
0082   const TrackerGeometry* pDD = &es.getData(tkGeomToken_);
0083   edm::LogInfo("SiPixelLorentzAngle (old)")
0084       << " There are " << pDD->detUnits().size() << " detectors (old)" << std::endl;
0085 
0086   for (const auto& it : pDD->detUnits()) {
0087     if (dynamic_cast<PixelGeomDetUnit const*>(it) != nullptr) {
0088       DetId detid = it->geographicalId();
0089       const DetId detidc = it->geographicalId();
0090 
0091       // fill bpix values for LA
0092       if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
0093         edm::LogPrint("SiPixelLorentzAngleDB")
0094             << " pixel barrel:"
0095             << "  layer=" << tTopo->pxbLayer(detidc.rawId()) << "  ladder=" << tTopo->pxbLadder(detidc.rawId())
0096             << "  module=" << tTopo->pxbModule(detidc.rawId()) << "  rawId=" << detidc.rawId() << endl;
0097 
0098         if (!useFile_) {
0099           //first individuals are put
0100           for (Parameters::iterator it = ModuleParameters_.begin(); it != ModuleParameters_.end(); ++it) {
0101             if (it->getParameter<unsigned int>("rawid") == detidc.rawId()) {
0102               float lorentzangle = (float)it->getParameter<double>("angle");
0103               LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0104               edm::LogPrint("SiPixelLorentzAngleDB")
0105                   << " individual value=" << lorentzangle << " put into rawid=" << detid.rawId() << endl;
0106             }
0107           }
0108 
0109           //modules already put are automatically skipped
0110           for (Parameters::iterator it = BPixParameters_.begin(); it != BPixParameters_.end(); ++it) {
0111             if (it->getParameter<unsigned int>("module") == tTopo->pxbModule(detidc.rawId()) &&
0112                 it->getParameter<unsigned int>("layer") == tTopo->pxbLayer(detidc.rawId())) {
0113               float lorentzangle = (float)it->getParameter<double>("angle");
0114               LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0115             }
0116           }
0117 
0118         } else {
0119           edm::LogError("SiPixelLorentzAngleDB")
0120               << "[SiPixelLorentzAngleDB::analyze] method for reading file not implemented yet" << std::endl;
0121         }
0122 
0123         // fill fpix values for LA
0124       } else if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
0125         edm::LogPrint("SiPixelLorentzAngleDB")
0126             << " pixel endcap:"
0127             << "  side=" << tTopo->pxfSide(detidc.rawId()) << "  disk=" << tTopo->pxfDisk(detidc.rawId())
0128             << "  blade=" << tTopo->pxfBlade(detidc.rawId()) << "  panel=" << tTopo->pxfPanel(detidc.rawId())
0129             << "  module=" << tTopo->pxfModule(detidc.rawId()) << "  rawId=" << detidc.rawId() << endl;
0130 
0131         //first individuals are put
0132         for (Parameters::iterator it = ModuleParameters_.begin(); it != ModuleParameters_.end(); ++it) {
0133           if (it->getParameter<unsigned int>("rawid") == detidc.rawId()) {
0134             float lorentzangle = (float)it->getParameter<double>("angle");
0135             LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0136             edm::LogPrint("SiPixelLorentzAngleDB")
0137                 << " individual value=" << lorentzangle << " put into rawid=" << detid.rawId() << endl;
0138           }
0139         }
0140 
0141         //modules already put are automatically skipped
0142         for (Parameters::iterator it = FPixParameters_.begin(); it != FPixParameters_.end(); ++it) {
0143           if (it->getParameter<unsigned int>("side") == tTopo->pxfSide(detidc.rawId()) &&
0144               it->getParameter<unsigned int>("disk") == tTopo->pxfDisk(detidc.rawId()) &&
0145               it->getParameter<unsigned int>("HVgroup") ==
0146                   HVgroup(tTopo->pxfPanel(detidc.rawId()), tTopo->pxfModule(detidc.rawId()))) {
0147             float lorentzangle = (float)it->getParameter<double>("angle");
0148             LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0149           }
0150         }
0151 
0152       } else {
0153         edm::LogError("SiPixelLorentzAngleDB")
0154             << "[SiPixelLorentzAngleDB::analyze] detid is Pixel but neither bpix nor fpix" << std::endl;
0155       }
0156     }
0157   }
0158 
0159   edm::Service<cond::service::PoolDBOutputService> mydbservice;
0160   if (mydbservice.isAvailable()) {
0161     try {
0162       if (mydbservice->isNewTagRequest(recordName_)) {
0163         mydbservice->createOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->beginOfTime(), recordName_);
0164       } else {
0165         mydbservice->appendOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->currentTime(), recordName_);
0166       }
0167     } catch (const cond::Exception& er) {
0168       edm::LogError("SiPixelLorentzAngleDB") << er.what() << std::endl;
0169     } catch (const std::exception& er) {
0170       edm::LogError("SiPixelLorentzAngleDB") << "caught std::exception " << er.what() << std::endl;
0171     } catch (...) {
0172       edm::LogError("SiPixelLorentzAngleDB") << "Funny error" << std::endl;
0173     }
0174   } else {
0175     edm::LogError("SiPixelLorentzAngleDB") << "Service is unavailable" << std::endl;
0176   }
0177 }
0178 
0179 unsigned int SiPixelLorentzAngleDB::HVgroup(unsigned int panel, unsigned int module) {
0180   if (1 == panel && (1 == module || 2 == module)) {
0181     return 1;
0182   } else if (1 == panel && (3 == module || 4 == module)) {
0183     return 2;
0184   } else if (2 == panel && 1 == module) {
0185     return 1;
0186   } else if (2 == panel && (2 == module || 3 == module)) {
0187     return 2;
0188   } else {
0189     edm::LogPrint("SiPixelLorentzAngleDB") << " *** error *** in SiPixelLorentzAngleDB::HVgroup(...), panel = " << panel
0190                                            << ", module = " << module << endl;
0191     return 0;
0192   }
0193 }
0194 DEFINE_FWK_MODULE(SiPixelLorentzAngleDB);