Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <memory>
0002 #include <string>
0003 #include <iostream>
0004 #include <fstream>
0005 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0006 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
0007 #include "DataFormats/DetId/interface/DetId.h"
0008 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementError.h"
0009 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
0010 #include "DataFormats/GeometrySurface/interface/GloballyPositioned.h"
0011 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.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/CommonTopologies/interface/PixelGeomDetUnit.h"
0023 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0024 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0025 
0026 class SiPixelLorentzAngleDBLoader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0027 public:
0028   explicit SiPixelLorentzAngleDBLoader(const edm::ParameterSet& conf);
0029 
0030   ~SiPixelLorentzAngleDBLoader() override = default;
0031   void analyze(const edm::Event& e, const edm::EventSetup& c) override;
0032 
0033 private:
0034   typedef std::vector<edm::ParameterSet> Parameters;
0035   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
0036   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tkTopoToken_;
0037   const std::string recordName_;
0038   const bool useFile_;
0039   const std::string fileName_;
0040   const Parameters BPixParameters_;
0041   const Parameters FPixParameters_;
0042   const Parameters ModuleParameters_;
0043   float bPixLorentzAnglePerTesla_;
0044   float fPixLorentzAnglePerTesla_;
0045 
0046   int HVgroup(int panel, int module);
0047 
0048   std::vector<std::pair<uint32_t, float> > detid_la;
0049 };
0050 
0051 using namespace std;
0052 using namespace edm;
0053 
0054 // Constructor
0055 SiPixelLorentzAngleDBLoader::SiPixelLorentzAngleDBLoader(edm::ParameterSet const& conf)
0056     : tkGeomToken_(esConsumes()),
0057       tkTopoToken_(esConsumes()),
0058       recordName_(conf.getUntrackedParameter<std::string>("record", "SiPixelLorentzAngleRcd")),
0059       useFile_(conf.getParameter<bool>("useFile")),
0060       fileName_(conf.getParameter<string>("fileName")),
0061       BPixParameters_(conf.getUntrackedParameter<Parameters>("BPixParameters")),
0062       FPixParameters_(conf.getUntrackedParameter<Parameters>("FPixParameters")),
0063       ModuleParameters_(conf.getUntrackedParameter<Parameters>("ModuleParameters")) {
0064   bPixLorentzAnglePerTesla_ =
0065       static_cast<float>(conf.getUntrackedParameter<double>("bPixLorentzAnglePerTesla", -9999.));
0066   fPixLorentzAnglePerTesla_ =
0067       static_cast<float>(conf.getUntrackedParameter<double>("fPixLorentzAnglePerTesla", -9999.));
0068   usesResource(cond::service::PoolDBOutputService::kSharedResource);
0069 }
0070 
0071 void SiPixelLorentzAngleDBLoader::analyze(const edm::Event& e, const edm::EventSetup& es) {
0072   static constexpr int nModules_ = 4;
0073   SiPixelLorentzAngle LorentzAngle;
0074 
0075   // Retrieve tracker geometry from geometry
0076   const TrackerGeometry* pDD = &es.getData(tkGeomToken_);
0077   // Retrieve tracker topology from geometry
0078   const TrackerTopology* tTopo = &es.getData(tkTopoToken_);
0079 
0080   for (auto& unit : pDD->detUnits()) {
0081     if (auto pixelUnit = dynamic_cast<PixelGeomDetUnit const*>(unit)) {
0082       const DetId detid = pixelUnit->geographicalId();
0083       auto rawId = detid.rawId();
0084       int found = 0;
0085       int side = tTopo->side(detid);  // 1:-z 2:+z for fpix, for bpix gives 0
0086 
0087       // fill bpix values for LA
0088       if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
0089         int layer = tTopo->pxbLayer(detid);
0090         // Barrel ladder id 1-20,32,44.
0091         int ladder = tTopo->pxbLadder(detid);
0092         // Barrel Z-index=1,8
0093         int module = tTopo->pxbModule(detid);
0094         if (module < nModules_ + 1) {
0095           side = 1;
0096         } else {
0097           side = 2;
0098         }
0099 
0100         LogPrint("SiPixelLorentzAngleDBLoader") << " pixel barrel:"
0101                                                 << " layer=" << layer << " ladder=" << ladder << " module=" << module
0102                                                 << "  rawId=" << rawId << " side=" << side;
0103 
0104         // use a commmon value (e.g. for MC)
0105         if (bPixLorentzAnglePerTesla_ != -9999.) {  // use common value for all
0106           LogPrint("SiPixelLorentzAngleDBLoader")
0107               << " LA=" << bPixLorentzAnglePerTesla_ << " common for all bpix" << endl;
0108           if (!LorentzAngle.putLorentzAngle(detid.rawId(), bPixLorentzAnglePerTesla_)) {
0109             LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
0110           }
0111           // use an external file
0112         } else if (useFile_) {
0113           LogPrint("SiPixelLorentzAngleDBLoader") << "method for reading file not implemented yet";
0114           // use config file
0115         } else {
0116           // first individuals are put
0117           for (auto& moduleParam : ModuleParameters_) {
0118             if (moduleParam.getParameter<unsigned int>("rawid") == detid.rawId()) {
0119               float lorentzangle = static_cast<float>(moduleParam.getParameter<double>("angle"));
0120               if (!found) {
0121                 LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0122                 LogPrint("SiPixelLorentzAngleDBLoader")
0123                     << "   >> LA=" << lorentzangle << " individual value " << detid.rawId();
0124                 found = 1;
0125               } else {
0126                 LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
0127               }
0128             }
0129           }  // end on loop for ModuleParameters_
0130 
0131           //modules already put are automatically skipped
0132           for (auto& bpixParam : BPixParameters_) {
0133             if (bpixParam.exists("layer")) {
0134               if (bpixParam.getParameter<int>("layer") != layer)
0135                 continue;
0136               if (bpixParam.exists("ladder"))
0137                 if (bpixParam.getParameter<int>("ladder") != ladder)
0138                   continue;
0139               if (bpixParam.exists("module"))
0140                 if (bpixParam.getParameter<int>("module") != module)
0141                   continue;
0142               if (bpixParam.exists("side"))
0143                 if (bpixParam.getParameter<int>("side") != side)
0144                   continue;
0145               if (!found) {
0146                 float lorentzangle = static_cast<float>(bpixParam.getParameter<double>("angle"));
0147                 LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0148                 LogPrint("SiPixelLorentzAngleDBLoader") << "   >> LA=" << lorentzangle;
0149                 found = 2;
0150               } else if (found == 1) {
0151                 LogPrint("SiPixelLorentzAngleDBLoader") << "The detid already given in ModuleParameters, skipping ...";
0152               } else
0153                 LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
0154             }
0155           }
0156         }  // condition to read from config
0157 
0158         // fill fpix values for LA (for phase2 fpix & epix)
0159       } else if (detid.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
0160         // Convert to online
0161         PixelEndcapName pen(detid, tTopo, true);
0162         int disk = pen.diskName();
0163         int blade = pen.bladeName();
0164         int panel = pen.pannelName();
0165         int ring = pen.ringName();
0166 
0167         LogPrint("SiPixelLorentzAngleDBLoader") << " pixel endcap:"
0168                                                 << " side=" << side << " disk=" << disk << " blade=" << blade
0169                                                 << " pannel=" << panel << " ring=" << ring << "  rawId=" << rawId;
0170 
0171         // use a commmon value (e.g. for MC)
0172         if (fPixLorentzAnglePerTesla_ != -9999.) {  // use common value for all
0173           LogPrint("SiPixelLorentzAngleDBLoader") << " LA =" << fPixLorentzAnglePerTesla_ << " common for all FPix";
0174           if (!LorentzAngle.putLorentzAngle(detid.rawId(), fPixLorentzAnglePerTesla_)) {
0175             LogError("SiPixelLorentzAngleDBLoader") << "detid already exists";
0176           }
0177 
0178         } else if (useFile_) {
0179           LogPrint("SiPixelLorentzAngleDBLoader") << "method for reading file not implemented yet";
0180 
0181         } else {
0182           //first individuals are put
0183           for (auto& parameter : ModuleParameters_) {
0184             if (parameter.getParameter<unsigned int>("rawid") == detid.rawId()) {
0185               float lorentzangle = static_cast<float>(parameter.getParameter<double>("angle"));
0186               if (!found) {
0187                 LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0188                 LogPrint("SiPixelLorentzAngleDBLoader")
0189                     << " LA=" << lorentzangle << " individual value " << detid.rawId();
0190                 found = 1;
0191               } else
0192                 LogError("SiPixelLorentzAngleDBLoader") << "ERROR!: detid already exists";
0193             }
0194           }  // end loop on ModuleParameters_
0195 
0196           // modules already put are automatically skipped
0197           for (auto& fpixParam : FPixParameters_) {
0198             if (fpixParam.exists("side"))
0199               if (fpixParam.getParameter<int>("side") != side)
0200                 continue;
0201             if (fpixParam.exists("disk"))
0202               if (fpixParam.getParameter<int>("disk") != disk)
0203                 continue;
0204             if (fpixParam.exists("ring"))
0205               if (fpixParam.getParameter<int>("ring") != ring)
0206                 continue;
0207             if (fpixParam.exists("blade"))
0208               if (fpixParam.getParameter<int>("blade") != blade)
0209                 continue;
0210             if (fpixParam.exists("panel"))
0211               if (fpixParam.getParameter<int>("panel") != panel)
0212                 continue;
0213             if (fpixParam.exists("HVgroup"))
0214               if (fpixParam.getParameter<int>("HVgroup") != HVgroup(panel, ring))
0215                 continue;
0216             if (!found) {
0217               float lorentzangle = static_cast<float>(fpixParam.getParameter<double>("angle"));
0218               LorentzAngle.putLorentzAngle(detid.rawId(), lorentzangle);
0219               LogPrint("SiPixelLorentzAngleDBLoader") << "   >> LA=" << lorentzangle;
0220               found = 2;
0221             } else if (found == 1) {
0222               LogPrint("SiPixelLorentzAngleDBLoader") << "The detid already given in ModuleParameters, skipping ...";
0223             } else
0224               LogError("SiPixelLorentzAngleDBLoader") << " ERROR!: detid already exists";
0225           }  // end loop on FPixParameters_
0226         }    // condition to read from config
0227       }      // end on being barrel or endcap
0228     }
0229   }
0230 
0231   edm::Service<cond::service::PoolDBOutputService> mydbservice;
0232   if (mydbservice.isAvailable()) {
0233     try {
0234       if (mydbservice->isNewTagRequest(recordName_)) {
0235         mydbservice->createOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->beginOfTime(), recordName_);
0236       } else {
0237         mydbservice->appendOneIOV<SiPixelLorentzAngle>(LorentzAngle, mydbservice->currentTime(), recordName_);
0238       }
0239     } catch (const cond::Exception& er) {
0240       LogPrint("SiPixelLorentzAngleDBLoader") << "SiPixelLorentzAngleDBLoader" << er.what();
0241     } catch (const std::exception& er) {
0242       LogPrint("SiPixelLorentzAngleDBLoader") << "SiPixelLorentzAngleDBLoader"
0243                                               << "caught std::exception " << er.what();
0244     }
0245   } else {
0246     LogPrint("SiPixelLorentzAngleDBLoader") << "Service is unavailable";
0247   }
0248 }
0249 
0250 int SiPixelLorentzAngleDBLoader::HVgroup(int panel, int module) {
0251   if (1 == panel && (1 == module || 2 == module)) {
0252     return 1;
0253   } else if (1 == panel && (3 == module || 4 == module)) {
0254     return 2;
0255   } else if (2 == panel && 1 == module) {
0256     return 1;
0257   } else if (2 == panel && (2 == module || 3 == module)) {
0258     return 2;
0259   } else {
0260     LogError("SiPixelLorentzAngleDBLoader")
0261         << " *** error *** in SiPixelLorentzAngleDBLoader::HVgroup(...), panel = " << panel << ", module = " << module
0262         << endl;
0263     return 0;
0264   }
0265 }
0266 
0267 //define this as a plug-in
0268 DEFINE_FWK_MODULE(SiPixelLorentzAngleDBLoader);