Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:02:25

0001 // -*- C++ -*-
0002 //
0003 // Package:    CondFormats/PCLConfig
0004 // Class:      AlignPCLThresholdsWriter
0005 //
0006 /**\class AlignPCLThresholdsWriter AlignPCLThresholdsWriter.cc CondFormats/PCLConfig/plugins/AlignPCLThresholdsWriter.cc
0007 
0008  Description: class to build the SiPixelAli PCL thresholds
0009 
0010 */
0011 //
0012 // Original Author:  Marco Musich
0013 //         Created:  Wed, 22 Feb 2017 12:04:36 GMT
0014 //
0015 //
0016 
0017 // system include files
0018 #include <array>
0019 #include <memory>
0020 
0021 // user include files
0022 #include "FWCore/Framework/interface/Frameworkfwd.h"
0023 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0024 
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027 
0028 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0029 #include "CondFormats/PCLConfig/interface/AlignPCLThresholds.h"
0030 #include "CondFormats/PCLConfig/interface/AlignPCLThresholdsHG.h"
0031 #include "FWCore/ServiceRegistry/interface/Service.h"
0032 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0033 
0034 //
0035 // class declaration
0036 //
0037 
0038 namespace DOFs {
0039   enum dof { X, Y, Z, thetaX, thetaY, thetaZ, extraDOF };
0040 }
0041 
0042 template <typename T>
0043 class AlignPCLThresholdsWriter : public edm::one::EDAnalyzer<> {
0044 public:
0045   explicit AlignPCLThresholdsWriter(const edm::ParameterSet&);
0046   ~AlignPCLThresholdsWriter() override = default;
0047 
0048   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0049 
0050 private:
0051   void analyze(const edm::Event&, const edm::EventSetup&) override;
0052   DOFs::dof mapOntoEnum(std::string coord);
0053 
0054   void writePayload(T& myThresholds);
0055   void storeHGthresholds(AlignPCLThresholdsHG& myThresholds, const std::vector<std::string>& alignables);
0056 
0057   // ----------member data ---------------------------
0058   const std::string m_record;
0059   const unsigned int m_minNrecords;
0060   const std::vector<edm::ParameterSet> m_parameters;
0061 };
0062 
0063 //
0064 // constructors and destructor
0065 //
0066 template <typename T>
0067 AlignPCLThresholdsWriter<T>::AlignPCLThresholdsWriter(const edm::ParameterSet& iConfig)
0068     : m_record(iConfig.getParameter<std::string>("record")),
0069       m_minNrecords(iConfig.getParameter<unsigned int>("minNRecords")),
0070       m_parameters(iConfig.getParameter<std::vector<edm::ParameterSet> >("thresholds")) {}
0071 
0072 //
0073 // member functions
0074 //
0075 
0076 // ------------ method called for each event  ------------
0077 template <typename T>
0078 void AlignPCLThresholdsWriter<T>::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0079   // detect if new payload is used
0080   bool newClass = false;
0081   for (auto& thePSet : m_parameters) {
0082     if (thePSet.exists("fractionCut")) {
0083       newClass = true;
0084       break;
0085     }
0086   }
0087 
0088   T myThresholds{};
0089   if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
0090     if (newClass) {
0091       this->writePayload(myThresholds);
0092     } else {
0093       throw cms::Exception("AlignPCLThresholdsWriter") << "mismatched configuration";
0094     }
0095   } else {
0096     if (!newClass) {
0097       this->writePayload(myThresholds);
0098     } else {
0099       throw cms::Exception("AlignPCLThresholdsWriter") << "mismatched configuration";
0100     }
0101   }
0102 }
0103 
0104 template <typename T>
0105 DOFs::dof AlignPCLThresholdsWriter<T>::mapOntoEnum(std::string coord) {
0106   if (coord == "X") {
0107     return DOFs::X;
0108   } else if (coord == "Y") {
0109     return DOFs::Y;
0110   } else if (coord == "Z") {
0111     return DOFs::Z;
0112   } else if (coord == "thetaX") {
0113     return DOFs::thetaX;
0114   } else if (coord == "thetaY") {
0115     return DOFs::thetaY;
0116   } else if (coord == "thetaZ") {
0117     return DOFs::thetaZ;
0118   } else {
0119     return DOFs::extraDOF;
0120   }
0121 }
0122 
0123 // ------------ templated method to write the payload  ------------
0124 template <typename T>
0125 void AlignPCLThresholdsWriter<T>::writePayload(T& myThresholds) {
0126   using namespace edm;
0127 
0128   edm::LogInfo("AlignPCLThresholdsWriter") << "Size of AlignPCLThresholds object " << myThresholds.size() << std::endl;
0129 
0130   // loop on the PSet and insert the conditions
0131   std::array<std::string, 6> mandatories = {{"X", "Y", "Z", "thetaX", "thetaY", "thetaZ"}};
0132   std::vector<std::string> alignables;
0133 
0134   // fill the list of alignables
0135   for (auto& thePSet : m_parameters) {
0136     const std::string alignableId(thePSet.getParameter<std::string>("alignableId"));
0137     // only if it is not yet in the list
0138     if (std::find(alignables.begin(), alignables.end(), alignableId) == alignables.end()) {
0139       alignables.push_back(alignableId);
0140     }
0141   }
0142 
0143   for (auto& alignable : alignables) {
0144     AlignPCLThreshold::coordThresholds my_X;
0145     AlignPCLThreshold::coordThresholds my_Y;
0146     AlignPCLThreshold::coordThresholds my_Z;
0147     AlignPCLThreshold::coordThresholds my_tX;
0148     AlignPCLThreshold::coordThresholds my_tY;
0149     AlignPCLThreshold::coordThresholds my_tZ;
0150 
0151     std::vector<std::string> presentDOF;
0152 
0153     // extra degrees of freedom
0154     std::vector<AlignPCLThreshold::coordThresholds> extraDOFs = std::vector<AlignPCLThreshold::coordThresholds>();
0155 
0156     for (auto& thePSet : m_parameters) {
0157       const std::string alignableId(thePSet.getParameter<std::string>("alignableId"));
0158       const std::string DOF(thePSet.getParameter<std::string>("DOF"));
0159 
0160       const double cutoff(thePSet.getParameter<double>("cut"));
0161       const double sigCut(thePSet.getParameter<double>("sigCut"));
0162       const double maxMoveCut(thePSet.getParameter<double>("maxMoveCut"));
0163       const double maxErrorCut(thePSet.getParameter<double>("maxErrorCut"));
0164 
0165       if (alignableId == alignable) {
0166         presentDOF.push_back(DOF);
0167         // create the objects
0168 
0169         switch (mapOntoEnum(DOF)) {
0170           case DOFs::X:
0171             my_X.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
0172             break;
0173           case DOFs::Y:
0174             my_Y.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
0175             break;
0176           case DOFs::Z:
0177             my_Z.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
0178             break;
0179           case DOFs::thetaX:
0180             my_tX.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
0181             break;
0182           case DOFs::thetaY:
0183             my_tY.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
0184             break;
0185           case DOFs::thetaZ:
0186             my_tZ.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
0187             break;
0188           default:
0189             edm::LogInfo("AlignPCLThresholdsWriter")
0190                 << "Appending Extra degree of freeedom: " << DOF << " " << mapOntoEnum(DOF) << std::endl;
0191             AlignPCLThreshold::coordThresholds ExtraDOF;
0192             ExtraDOF.setThresholds(cutoff, sigCut, maxErrorCut, maxMoveCut, DOF);
0193             extraDOFs.push_back(ExtraDOF);
0194         }
0195 
0196         AlignPCLThreshold a(my_X, my_tX, my_Y, my_tY, my_Z, my_tZ, extraDOFs);
0197         myThresholds.setAlignPCLThreshold(alignableId, a);
0198 
0199       }  // if alignable is found in the PSet
0200     }    // loop on the PSets
0201 
0202     // checks if all mandatories are present
0203     edm::LogInfo("AlignPCLThresholdsWriter")
0204         << "Size of AlignPCLThresholds object  " << myThresholds.size() << std::endl;
0205     for (auto& mandatory : mandatories) {
0206       if (std::find(presentDOF.begin(), presentDOF.end(), mandatory) == presentDOF.end()) {
0207         edm::LogWarning("AlignPCLThresholdsWriter")
0208             << "Configuration for DOF: " << mandatory << " for alignable " << alignable << "is not present \n"
0209             << "Will build object with defaults!" << std::endl;
0210       }
0211     }
0212 
0213   }  // ends loop on the alignable units
0214 
0215   // set the minimum number of records to be used in pede
0216   myThresholds.setNRecords(m_minNrecords);
0217   edm::LogInfo("AlignPCLThresholdsWriter") << "Content of AlignPCLThresholds " << std::endl;
0218 
0219   // additional thresholds for AlignPCLThresholdsHG
0220   if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
0221     storeHGthresholds(myThresholds, alignables);
0222   }
0223 
0224   // use built-in method in the CondFormat
0225   myThresholds.printAll();
0226 
0227   // Form the data here
0228   edm::Service<cond::service::PoolDBOutputService> poolDbService;
0229   if (poolDbService.isAvailable()) {
0230     cond::Time_t valid_time = poolDbService->currentTime();
0231     // this writes the payload to begin in current run defined in cfg
0232     poolDbService->writeOneIOV(myThresholds, valid_time, m_record);
0233   }
0234 }
0235 
0236 // ------------ method to store additional HG thresholds ------------
0237 template <typename T>
0238 void AlignPCLThresholdsWriter<T>::storeHGthresholds(AlignPCLThresholdsHG& myThresholds,
0239                                                     const std::vector<std::string>& alignables) {
0240   edm::LogInfo("AlignPCLThresholdsWriter")
0241       << "Found type AlignPCLThresholdsHG, additional thresholds are written" << std::endl;
0242 
0243   for (auto& alignable : alignables) {
0244     for (auto& thePSet : m_parameters) {
0245       const std::string alignableId(thePSet.getParameter<std::string>("alignableId"));
0246       const std::string DOF(thePSet.getParameter<std::string>("DOF"));
0247 
0248       // Get coordType from DOF
0249       AlignPCLThresholds::coordType type = static_cast<AlignPCLThresholds::coordType>(mapOntoEnum(DOF));
0250 
0251       if (alignableId == alignable) {
0252         if (thePSet.exists("fractionCut")) {
0253           const double fractionCut(thePSet.getParameter<double>("fractionCut"));
0254           myThresholds.setFractionCut(alignableId, type, fractionCut);
0255         }
0256       }
0257     }
0258   }
0259 }
0260 
0261 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0262 template <typename T>
0263 void AlignPCLThresholdsWriter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0264   edm::ParameterSetDescription desc;
0265   desc.setComment("Plugin to write payloads of type AlignPCLThresholds");
0266   desc.add<unsigned int>("minNRecords", 25000);
0267   edm::ParameterSetDescription desc_thresholds;
0268 
0269   desc_thresholds.add<std::string>("alignableId");
0270   desc_thresholds.add<std::string>("DOF");
0271   desc_thresholds.add<double>("cut");
0272   desc_thresholds.add<double>("sigCut");
0273   desc_thresholds.add<double>("maxMoveCut");
0274   desc_thresholds.add<double>("maxErrorCut");
0275   if constexpr (std::is_same_v<T, AlignPCLThresholdsHG>) {
0276     desc.add<std::string>("record", "AlignPCLThresholdsHGRcd");
0277     //optional thresholds from new payload version (not for all the alignables)
0278     desc_thresholds.addOptional<double>("fractionCut");
0279   } else {
0280     desc.add<std::string>("record", "AlignPCLThresholdsRcd");
0281   }
0282 
0283   std::vector<edm::ParameterSet> default_thresholds(1);
0284   desc.addVPSet("thresholds", desc_thresholds, default_thresholds);
0285   descriptions.addWithDefaultLabel(desc);
0286 }
0287 
0288 typedef AlignPCLThresholdsWriter<AlignPCLThresholds> AlignPCLThresholdsLGWriter;
0289 typedef AlignPCLThresholdsWriter<AlignPCLThresholdsHG> AlignPCLThresholdsHGWriter;
0290 
0291 DEFINE_FWK_MODULE(AlignPCLThresholdsLGWriter);
0292 DEFINE_FWK_MODULE(AlignPCLThresholdsHGWriter);