Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-10 23:04:48

0001 /*** Header file ***/
0002 #include "Alignment/MillePedeAlignmentAlgorithm/interface/MillePedeFileReader.h"
0003 #include "Alignment/CommonAlignment/interface/AlignableObjectId.h"
0004 
0005 /*** system includes ***/
0006 #include <cmath>  // include floating-point std::abs functions
0007 #include <fstream>
0008 
0009 /*** Alignment ***/
0010 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 
0013 //=============================================================================
0014 //===   PUBLIC METHOD IMPLEMENTATION                                        ===
0015 //=============================================================================
0016 
0017 MillePedeFileReader ::MillePedeFileReader(const edm::ParameterSet& config,
0018                                           const std::shared_ptr<const PedeLabelerBase>& pedeLabeler,
0019                                           const std::shared_ptr<const AlignPCLThresholdsHG>& theThresholds,
0020                                           const std::shared_ptr<const PixelTopologyMap>& pixelTopologyMap,
0021                                           const std::shared_ptr<const SiPixelQuality>& pixelQualityMap)
0022     : pedeLabeler_(pedeLabeler),
0023       theThresholds_(theThresholds),
0024       pixelTopologyMap_(pixelTopologyMap),
0025       ignoreInactiveAlignables_(config.getParameter<bool>("ignoreInactiveAlignables")),
0026       quality_(pixelQualityMap),
0027       dirName_(config.getParameter<std::string>("fileDir")),
0028       millePedeEndFile_(config.getParameter<std::string>("millePedeEndFile")),
0029       millePedeLogFile_(config.getParameter<std::string>("millePedeLogFile")),
0030       millePedeResFile_(config.getParameter<std::string>("millePedeResFile")),
0031       isHG_(config.getParameter<bool>("isHG")) {
0032   if (!dirName_.empty() && dirName_.find_last_of('/') != dirName_.size() - 1)
0033     dirName_ += '/';  // may need '/'
0034 }
0035 
0036 void MillePedeFileReader ::read() {
0037   if (isHG_) {
0038     initializeIndexHelper();
0039   }
0040   readMillePedeEndFile();
0041   readMillePedeLogFile();
0042   readMillePedeResultFile();
0043 }
0044 
0045 bool MillePedeFileReader ::storeAlignments() { return (updateDB_ && !vetoUpdateDB_); }
0046 
0047 //=============================================================================
0048 //===   PRIVATE METHOD IMPLEMENTATION                                       ===
0049 //=============================================================================
0050 void MillePedeFileReader ::readMillePedeEndFile() {
0051   std::ifstream endFile;
0052   endFile.open((dirName_ + millePedeEndFile_).c_str());
0053 
0054   if (endFile.is_open()) {
0055     edm::LogInfo("MillePedeFileReader") << "Reading millepede end-file";
0056     std::string line;
0057     getline(endFile, line);
0058     std::string trash;
0059     if (line.find("-1") != std::string::npos) {
0060       getline(endFile, line);
0061       exitMessage_ = line;
0062       std::istringstream iss(line);
0063       iss >> exitCode_ >> trash;
0064       edm::LogInfo("MillePedeFileReader")
0065           << " Pede exit code is: " << exitCode_ << " (" << exitMessage_ << ")" << std::endl;
0066     } else {
0067       exitMessage_ = line;
0068       std::istringstream iss(line);
0069       iss >> exitCode_ >> trash;
0070       edm::LogInfo("MillePedeFileReader")
0071           << " Pede exit code is: " << exitCode_ << " (" << exitMessage_ << ")" << std::endl;
0072     }
0073   } else {
0074     edm::LogError("MillePedeFileReader") << "Could not read millepede end-file.";
0075     exitMessage_ = "no exit code found";
0076   }
0077 }
0078 
0079 void MillePedeFileReader ::readMillePedeLogFile() {
0080   std::ifstream logFile;
0081   logFile.open((dirName_ + millePedeLogFile_).c_str());
0082 
0083   if (logFile.is_open()) {
0084     edm::LogInfo("MillePedeFileReader") << "Reading millepede log-file";
0085     std::string line;
0086 
0087     while (getline(logFile, line)) {
0088       std::string Nrec_string = "NREC =";
0089       std::string Binaries_string = "C_binary";
0090 
0091       if (line.find(Nrec_string) != std::string::npos) {
0092         std::istringstream iss(line);
0093         std::string trash;
0094         iss >> trash >> trash >> Nrec_;
0095 
0096         if (Nrec_ < theThresholds_->getNrecords()) {
0097           edm::LogInfo("MillePedeFileReader")
0098               << "Number of records used " << theThresholds_->getNrecords() << std::endl;
0099           updateDB_ = false;
0100         }
0101       }
0102 
0103       if (line.find(Binaries_string) != std::string::npos) {
0104         binariesAmount_ += 1;
0105       }
0106     }
0107   } else {
0108     edm::LogError("MillePedeFileReader") << "Could not read millepede log-file.";
0109 
0110     updateDB_ = false;
0111     Nrec_ = 0;
0112   }
0113 }
0114 
0115 void MillePedeFileReader ::readMillePedeResultFile() {
0116   // cutoffs by coordinate and by alignable
0117   std::map<std::string, std::array<float, 6> > cutoffs_;
0118   std::map<std::string, std::array<float, 6> > significances_;
0119   std::map<std::string, std::array<float, 6> > thresholds_;
0120   std::map<std::string, std::array<float, 6> > errors_;
0121   std::map<std::string, std::array<float, 6> > fractions_;
0122 
0123   std::map<std::string, std::array<int, 6> > countsAbove_;
0124   std::map<std::string, std::array<int, 6> > countsTotal_;
0125 
0126   AlignableObjectId alignableObjectId{AlignableObjectId::Geometry::General};
0127 
0128   std::vector<std::string> alignables_ = theThresholds_->getAlignableList();
0129   for (auto& ali : alignables_) {
0130     cutoffs_[ali] = theThresholds_->getCut(ali);
0131     significances_[ali] = theThresholds_->getSigCut(ali);
0132     thresholds_[ali] = theThresholds_->getMaxMoveCut(ali);
0133     errors_[ali] = theThresholds_->getMaxErrorCut(ali);
0134 
0135     if (theThresholds_->hasFloatMap(ali)) {
0136       fractions_[ali] = theThresholds_->getFractionCut(ali);
0137       countsAbove_[ali] = {{0, 0, 0, 0, 0, 0}};
0138       countsTotal_[ali] = {{0, 0, 0, 0, 0, 0}};
0139     }
0140   }
0141 
0142   updateDB_ = false;
0143   vetoUpdateDB_ = false;
0144   std::ifstream resFile;
0145   resFile.open((dirName_ + millePedeResFile_).c_str());
0146 
0147   if (resFile.is_open()) {
0148     edm::LogInfo("MillePedeFileReader") << "Reading millepede result-file";
0149 
0150     std::string line;
0151     getline(resFile, line);  // drop first line
0152 
0153     while (getline(resFile, line)) {
0154       std::istringstream iss(line);
0155 
0156       std::vector<std::string> tokens;
0157       std::string token;
0158       while (iss >> token) {
0159         tokens.push_back(token);
0160       }
0161 
0162       auto alignableLabel = std::stoul(tokens[0]);
0163       const auto alignable = pedeLabeler_->alignableFromLabel(alignableLabel);
0164       auto det = getHLS(alignable);
0165       // check if the modules associated to the alignable are active
0166       const bool active = isAlignableActive(alignable, quality_);
0167       int detIndex = static_cast<int>(det);
0168       auto alignableIndex = alignableLabel % 10 - 1;
0169       std::string detLabel = getStringFromHLS(det);
0170 
0171       if (!active) {
0172         edm::LogPrint("MillePedeFileReader") << "Alignable :" << detLabel << " is inactive";
0173       }
0174 
0175       if (tokens.size() > 4 /*3*/) {
0176         countsTotal_[detLabel][alignableIndex]++;  //Count aligned modules/ladders per structure
0177         const auto paramNum = pedeLabeler_->paramNumFromLabel(alignableLabel);
0178         align::StructureType type = alignable->alignableObjectId();
0179         align::ID id = alignable->id();
0180 
0181         double ObsMove = std::stof(tokens[3]) * multiplier_[alignableIndex];
0182         double ObsErr = std::stof(tokens[4]) * multiplier_[alignableIndex];
0183 
0184         auto coord = static_cast<AlignPCLThresholdsHG::coordType>(alignableIndex);
0185 
0186         if (det != PclHLS::NotInPCL) {
0187           if (type != align::TPBLadder && type != align::TPEPanel) {
0188             switch (coord) {
0189               case AlignPCLThresholdsHG::X:
0190                 Xobs_[detIndex] = ObsMove;
0191                 XobsErr_[detIndex] = ObsErr;
0192                 break;
0193               case AlignPCLThresholdsHG::Y:
0194                 Yobs_[detIndex] = ObsMove;
0195                 YobsErr_[detIndex] = ObsErr;
0196                 break;
0197               case AlignPCLThresholdsHG::Z:
0198                 Zobs_[detIndex] = ObsMove;
0199                 ZobsErr_[detIndex] = ObsErr;
0200                 break;
0201               case AlignPCLThresholdsHG::theta_X:
0202                 tXobs_[detIndex] = ObsMove;
0203                 tXobsErr_[detIndex] = ObsErr;
0204                 break;
0205               case AlignPCLThresholdsHG::theta_Y:
0206                 tYobs_[detIndex] = ObsMove;
0207                 tYobsErr_[detIndex] = ObsErr;
0208                 break;
0209               case AlignPCLThresholdsHG::theta_Z:
0210                 tZobs_[detIndex] = ObsMove;
0211                 tZobsErr_[detIndex] = ObsErr;
0212                 break;
0213               default:
0214                 edm::LogError("MillePedeFileReader") << "Currently not able to handle DOF " << coord << std::endl;
0215                 break;
0216             }
0217           } else {
0218             auto hgIndex = getIndexForHG(id, det);
0219             switch (coord) {
0220               case AlignPCLThresholdsHG::X:
0221                 Xobs_HG_[hgIndex - 1] = ObsMove;
0222                 XobsErr_HG_[hgIndex - 1] = ObsErr;
0223                 break;
0224               case AlignPCLThresholdsHG::Y:
0225                 Yobs_HG_[hgIndex - 1] = ObsMove;
0226                 YobsErr_HG_[hgIndex - 1] = ObsErr;
0227                 break;
0228               case AlignPCLThresholdsHG::Z:
0229                 Zobs_HG_[hgIndex - 1] = ObsMove;
0230                 ZobsErr_HG_[hgIndex - 1] = ObsErr;
0231                 break;
0232               case AlignPCLThresholdsHG::theta_X:
0233                 tXobs_HG_[hgIndex - 1] = ObsMove;
0234                 tXobsErr_HG_[hgIndex - 1] = ObsErr;
0235                 break;
0236               case AlignPCLThresholdsHG::theta_Y:
0237                 tYobs_HG_[hgIndex - 1] = ObsMove;
0238                 tYobsErr_HG_[hgIndex - 1] = ObsErr;
0239                 break;
0240               case AlignPCLThresholdsHG::theta_Z:
0241                 tZobs_HG_[hgIndex - 1] = ObsMove;
0242                 tZobsErr_HG_[hgIndex - 1] = ObsErr;
0243                 break;
0244               default:
0245                 edm::LogError("MillePedeFileReader") << "Currently not able to handle DOF " << coord << std::endl;
0246                 break;
0247             }
0248           }
0249 
0250         } else {
0251           edm::LogError("MillePedeFileReader")
0252               << "Currently not able to handle coordinate: " << coord << " (" << paramNum << ")  "
0253               << Form(" %s with ID %d (subdet %d)", alignableObjectId.idToString(type), id, DetId(id).subdetId())
0254               << std::endl;
0255           continue;
0256         }
0257 
0258         edm::LogVerbatim("MillePedeFileReader")
0259             << " alignableLabel: " << alignableLabel << " with alignableIndex " << alignableIndex << " detIndex "
0260             << detIndex << "\n"
0261             << " i.e. detLabel: " << detLabel << " (" << coord << ")\n"
0262             << " has movement: " << ObsMove << " +/- " << ObsErr << "\n"
0263             << " cutoff (cutoffs_[" << detLabel << "][" << coord << "]): " << cutoffs_[detLabel][alignableIndex] << "\n"
0264             << " significance (significances_[" << detLabel << "][" << coord
0265             << "]): " << significances_[detLabel][alignableIndex] << "\n"
0266             << " error thresolds (errors_[" << detLabel << "][" << coord << "]): " << errors_[detLabel][alignableIndex]
0267             << "\n"
0268             << " max movement (thresholds_[" << detLabel << "][" << coord
0269             << "]): " << thresholds_[detLabel][alignableIndex] << "\n"
0270             << " fraction (fractions_[" << detLabel << "][" << coord << "]): " << fractions_[detLabel][alignableIndex]
0271             << "\n"
0272             << "=============" << std::endl;
0273 
0274         if (std::abs(ObsMove) > thresholds_[detLabel][alignableIndex]) {
0275           if (active || ignoreInactiveAlignables_) {
0276             edm::LogWarning("MillePedeFileReader")
0277                 << "Aborting payload creation."
0278                 << " Exceeding maximum thresholds for movement: " << std::abs(ObsMove) << " for " << detLabel << " ("
0279                 << coord << ")";
0280             updateBits_.set(0);
0281             vetoUpdateDB_ = true;
0282             continue;
0283           } else {
0284             edm::LogInfo("MillePedeFileReader")
0285                 << " Exceeding maximum thresholds for movement: " << std::abs(ObsMove) << " for " << detLabel << " ("
0286                 << coord << ") but continuing as the alignable is inactive!";
0287           }
0288 
0289         } else if (std::abs(ObsMove) > cutoffs_[detLabel][alignableIndex]) {
0290           updateBits_.set(1);
0291 
0292           if (std::abs(ObsErr) > errors_[detLabel][alignableIndex]) {
0293             edm::LogWarning("MillePedeFileReader") << "Aborting payload creation."
0294                                                    << " Exceeding maximum thresholds for error: " << std::abs(ObsErr)
0295                                                    << " for" << detLabel << "(" << coord << ")";
0296             updateBits_.set(2);
0297             vetoUpdateDB_ = true;
0298             continue;
0299           } else {
0300             if (std::abs(ObsMove / ObsErr) < significances_[detLabel][alignableIndex]) {
0301               updateBits_.set(3);
0302               continue;
0303             }
0304           }
0305           updateDB_ = true;
0306           if (!isHG_) {
0307             edm::LogInfo("MillePedeFileReader")
0308                 << "This correction: " << ObsMove << "+/-" << ObsErr << " for " << detLabel << "(" << coord
0309                 << ") will trigger a new Tracker Alignment payload!";
0310           }
0311           countsAbove_[detLabel][alignableIndex]++;
0312         }
0313       }
0314     }
0315   } else {
0316     edm::LogError("MillePedeFileReader") << "Could not read millepede result-file.";
0317 
0318     updateDB_ = false;
0319     Nrec_ = 0;
0320   }
0321 
0322   if (isHG_) {          // check fractionCut
0323     updateDB_ = false;  // reset booleans since fractionCut is considered for HG
0324     std::stringstream ss;
0325     for (auto& ali : alignables_) {
0326       ss << ali << std::endl;
0327       for (long unsigned int i = 0; i < countsTotal_[ali].size(); i++) {
0328         if (countsTotal_[ali][i] != 0 && fractions_[ali][i] != -1) {
0329           float fraction_ = countsAbove_[ali][i] / (1.0 * countsTotal_[ali][i]);
0330           ss << static_cast<AlignPCLThresholdsHG::coordType>(i) << ":   Fraction = " << fraction_
0331              << "   Fraction Threshold = " << fractions_[ali][i];
0332           if (fraction_ >= fractions_[ali][i]) {
0333             updateDB_ = true;
0334             ss << "   above fraction threshold" << std::endl;
0335             fractionExceeded_[ali][i] = true;
0336           } else {
0337             ss << std::endl;
0338             fractionExceeded_[ali][i] = false;
0339           }
0340         } else
0341           ss << "No entries available or no fraction thresholds defined" << std::endl;
0342       }
0343       ss << "===================" << std::endl;
0344     }
0345     if (updateDB_ && !vetoUpdateDB_) {
0346       ss << "Alignment will be updated" << std::endl;
0347     } else {
0348       ss << "Alignment will NOT be updated" << std::endl;
0349     }
0350     edm::LogWarning("MillePedeFileReader") << ss.str();
0351   }
0352 }
0353 
0354 bool MillePedeFileReader::isAlignableActive(const Alignable* alignable,
0355                                             const std::shared_ptr<const SiPixelQuality>& pixelQual) {
0356   std::vector<DetId> detIds;
0357 
0358   // Get the list of deep components (lowest daughters) of the Alignable
0359   const auto& deepComponents = alignable->deepComponents();
0360 
0361   // Iterate through the deep components to retrieve their DetIds
0362   for (const auto& component : deepComponents) {
0363     DetId detId = component->geomDetId();
0364     if (detId != DetId(0)) {
0365       detIds.push_back(detId);
0366     }
0367   }
0368 
0369   // Counter for bad modules
0370   int badModuleCount = 0;
0371   int totalDetIds = detIds.size();
0372 
0373   const auto& theDisabledModules = pixelQual->getBadComponentList();
0374   std::vector<SiPixelQuality::disabledModuleType> theDisabledModuleInAlignable;
0375   for (const auto& mod : theDisabledModules) {
0376     if (std::find(detIds.begin(), detIds.end(), mod.DetID) != detIds.end()) {
0377       theDisabledModuleInAlignable.push_back(mod);
0378     }
0379   }
0380 
0381   // nothing left to do
0382   if (theDisabledModuleInAlignable.empty()) {
0383     return true;
0384   }
0385 
0386   bool header{false};
0387 
0388   std::stringstream out;
0389   for (const auto& mod : theDisabledModuleInAlignable) {
0390     if (!header)
0391       out << " Alignable = " << getStringFromHLS(getHLS(alignable));
0392     header = true;
0393     bool isBad = pixelQual->IsModuleBad(mod.DetID);
0394     std::bitset<16> bad_rocs(mod.BadRocs);
0395     if (isBad || bad_rocs.all()) {
0396       badModuleCount++;
0397     }
0398     out << " " << mod.DetID << " (" << (isBad || bad_rocs.all()) << ") , " << bad_rocs;
0399   }
0400   out << std::endl;
0401 
0402   if (badModuleCount > 0) {
0403     out << " " << badModuleCount << " modules are bad out of " << totalDetIds << std::endl;
0404     edm::LogPrint("MillePedeFileReader") << out.str();
0405   }
0406 
0407   // Return false if at least half of the detIds are bad
0408   if (badModuleCount >= (totalDetIds + 1) / 2) {
0409     return false;
0410   }
0411 
0412   // If less than half are bad, return true
0413   return true;
0414 }
0415 
0416 MillePedeFileReader::PclHLS MillePedeFileReader ::getHLS(const Alignable* alignable) {
0417   if (!alignable)
0418     return PclHLS::NotInPCL;
0419 
0420   const auto& tns = pedeLabeler_->alignableTracker()->trackerNameSpace();
0421   const align::ID id = alignable->id();
0422 
0423   switch (alignable->alignableObjectId()) {
0424     case align::TPBHalfBarrel:
0425       switch (tns.tpb().halfBarrelNumber(id)) {
0426         case 1:
0427           return PclHLS::TPBHalfBarrelXminus;
0428         case 2:
0429           return PclHLS::TPBHalfBarrelXplus;
0430         default:
0431           throw cms::Exception("LogicError")
0432               << "@SUB=MillePedeFileReader::getHLS\n"
0433               << "Found a pixel half-barrel number that should not exist: " << tns.tpb().halfBarrelNumber(id);
0434       }
0435     case align::TPEHalfCylinder:
0436       switch (tns.tpe().endcapNumber(id)) {
0437         case 1:
0438           switch (tns.tpe().halfCylinderNumber(id)) {
0439             case 1:
0440               return PclHLS::TPEHalfCylinderXminusZminus;
0441             case 2:
0442               return PclHLS::TPEHalfCylinderXplusZminus;
0443             default:
0444               throw cms::Exception("LogicError")
0445                   << "@SUB=MillePedeFileReader::getHLS\n"
0446                   << "Found a pixel half-cylinder number that should not exist: " << tns.tpe().halfCylinderNumber(id);
0447           }
0448         case 2:
0449           switch (tns.tpe().halfCylinderNumber(id)) {
0450             case 1:
0451               return PclHLS::TPEHalfCylinderXminusZplus;
0452             case 2:
0453               return PclHLS::TPEHalfCylinderXplusZplus;
0454             default:
0455               throw cms::Exception("LogicError")
0456                   << "@SUB=MillePedeFileReader::getHLS\n"
0457                   << "Found a pixel half-cylinder number that should not exist: " << tns.tpe().halfCylinderNumber(id);
0458           }
0459         default:
0460           throw cms::Exception("LogicError")
0461               << "@SUB=MillePedeFileReader::getHLS\n"
0462               << "Found a pixel endcap number that should not exist: " << tns.tpe().endcapNumber(id);
0463       }
0464     case align::TPBLadder:
0465       switch (tns.tpb().layerNumber(id)) {
0466         case 1:
0467           return PclHLS::TPBLadderLayer1;
0468         case 2:
0469           return PclHLS::TPBLadderLayer2;
0470         case 3:
0471           return PclHLS::TPBLadderLayer3;
0472         case 4:
0473           return PclHLS::TPBLadderLayer4;
0474         default:
0475           throw cms::Exception("LogicError")
0476               << "@SUB=MillePedeFileReader::getHLS\n"
0477               << "Found a pixel layer number that should not exist: " << tns.tpb().layerNumber(id);
0478       }
0479     case align::TPEPanel:
0480       switch (static_cast<signed int>((tns.tpe().endcapNumber(id) == 1) ? -1 * tns.tpe().halfDiskNumber(id)
0481                                                                         : tns.tpe().halfDiskNumber(id))) {
0482         case -3:
0483           return PclHLS::TPEPanelDiskM3;
0484         case -2:
0485           return PclHLS::TPEPanelDiskM2;
0486         case -1:
0487           return PclHLS::TPEPanelDiskM1;
0488         case 3:
0489           return PclHLS::TPEPanelDisk3;
0490         case 2:
0491           return PclHLS::TPEPanelDisk2;
0492         case 1:
0493           return PclHLS::TPEPanelDisk1;
0494         default:
0495           throw cms::Exception("LogicError")
0496               << "@SUB=MillePedeFileReader::getHLS\n"
0497               << "Found a pixel disk number that should not exist: "
0498               << static_cast<signed int>((tns.tpe().endcapNumber(id) == 1) ? -1 * tns.tpe().halfDiskNumber(id)
0499                                                                            : tns.tpe().halfDiskNumber(id));
0500       }
0501     default:
0502       return PclHLS::NotInPCL;
0503   }
0504 }
0505 
0506 std::string MillePedeFileReader::getStringFromHLS(MillePedeFileReader::PclHLS HLS) {
0507   switch (HLS) {
0508     case PclHLS::TPBHalfBarrelXminus:
0509       return "TPBHalfBarrelXminus";
0510     case PclHLS::TPBHalfBarrelXplus:
0511       return "TPBHalfBarrelXplus";
0512     case PclHLS::TPEHalfCylinderXminusZminus:
0513       return "TPEHalfCylinderXminusZminus";
0514     case PclHLS::TPEHalfCylinderXplusZminus:
0515       return "TPEHalfCylinderXplusZminus";
0516     case PclHLS::TPEHalfCylinderXminusZplus:
0517       return "TPEHalfCylinderXminusZplus";
0518     case PclHLS::TPEHalfCylinderXplusZplus:
0519       return "TPEHalfCylinderXplusZplus";
0520     case PclHLS::TPBLadderLayer1:
0521       return "TPBLadderLayer1";
0522     case PclHLS::TPBLadderLayer2:
0523       return "TPBLadderLayer2";
0524     case PclHLS::TPBLadderLayer3:
0525       return "TPBLadderLayer3";
0526     case PclHLS::TPBLadderLayer4:
0527       return "TPBLadderLayer4";
0528     case PclHLS::TPEPanelDisk1:
0529       return "TPEPanelDisk1";
0530     case PclHLS::TPEPanelDisk2:
0531       return "TPEPanelDisk2";
0532     case PclHLS::TPEPanelDisk3:
0533       return "TPEPanelDisk3";
0534     case PclHLS::TPEPanelDiskM1:
0535       return "TPEPanelDiskM1";
0536     case PclHLS::TPEPanelDiskM2:
0537       return "TPEPanelDiskM2";
0538     case PclHLS::TPEPanelDiskM3:
0539       return "TPEPanelDiskM3";
0540     default:
0541       //return "NotInPCL";
0542       throw cms::Exception("LogicError")
0543           << "@SUB=MillePedeFileReader::getStringFromHLS\n"
0544           << "Found an alignable structure not possible to map in the default AlignPCLThresholdsHG partitions";
0545   }
0546 }
0547 
0548 void MillePedeFileReader::initializeIndexHelper() {
0549   int currentSum = 0;
0550 
0551   indexHelper[PclHLS::TPBLadderLayer1] =
0552       std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXBLadders(1) / 2);
0553   currentSum += pixelTopologyMap_->getPXBLadders(1);
0554   indexHelper[PclHLS::TPBLadderLayer2] =
0555       std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXBLadders(2) / 2);
0556   currentSum += pixelTopologyMap_->getPXBLadders(2);
0557   indexHelper[PclHLS::TPBLadderLayer3] =
0558       std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXBLadders(3) / 2);
0559   currentSum += pixelTopologyMap_->getPXBLadders(3);
0560   indexHelper[PclHLS::TPBLadderLayer4] =
0561       std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXBLadders(4) / 2);
0562   currentSum += pixelTopologyMap_->getPXBLadders(4);
0563 
0564   indexHelper[PclHLS::TPEPanelDiskM3] = std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXFBlades(-3));
0565   currentSum += pixelTopologyMap_->getPXFBlades(-3) * 2;
0566   indexHelper[PclHLS::TPEPanelDiskM2] = std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXFBlades(-2));
0567   currentSum += pixelTopologyMap_->getPXFBlades(-2) * 2;
0568   indexHelper[PclHLS::TPEPanelDiskM1] = std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXFBlades(-1));
0569   currentSum += pixelTopologyMap_->getPXFBlades(-1) * 2;
0570 
0571   indexHelper[PclHLS::TPEPanelDisk1] = std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXFBlades(1));
0572   currentSum += pixelTopologyMap_->getPXFBlades(1) * 2;
0573   indexHelper[PclHLS::TPEPanelDisk2] = std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXFBlades(2));
0574   currentSum += pixelTopologyMap_->getPXFBlades(2) * 2;
0575   indexHelper[PclHLS::TPEPanelDisk3] = std::make_pair(currentSum, currentSum + pixelTopologyMap_->getPXFBlades(3));
0576 }
0577 
0578 int MillePedeFileReader::getIndexForHG(align::ID id, PclHLS HLS) {
0579   const auto& tns = pedeLabeler_->alignableTracker()->trackerNameSpace();
0580 
0581   switch (HLS) {
0582     case PclHLS::TPBLadderLayer1:
0583       return (tns.tpb().halfBarrelNumber(id) == 1) ? tns.tpb().ladderNumber(id) + indexHelper[HLS].first
0584                                                    : tns.tpb().ladderNumber(id) + indexHelper[HLS].second;
0585     case PclHLS::TPBLadderLayer2:
0586       return (tns.tpb().halfBarrelNumber(id) == 1) ? tns.tpb().ladderNumber(id) + indexHelper[HLS].first
0587                                                    : tns.tpb().ladderNumber(id) + indexHelper[HLS].second;
0588     case PclHLS::TPBLadderLayer3:
0589       return (tns.tpb().halfBarrelNumber(id) == 1) ? tns.tpb().ladderNumber(id) + indexHelper[HLS].first
0590                                                    : tns.tpb().ladderNumber(id) + indexHelper[HLS].second;
0591     case PclHLS::TPBLadderLayer4:
0592       return (tns.tpb().halfBarrelNumber(id) == 1) ? tns.tpb().ladderNumber(id) + indexHelper[HLS].first
0593                                                    : tns.tpb().ladderNumber(id) + indexHelper[HLS].second;
0594     case PclHLS::TPEPanelDisk1:
0595       return (tns.tpe().halfCylinderNumber(id) == 1)
0596                  ? (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].first
0597                  : (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].second;
0598     case PclHLS::TPEPanelDisk2:
0599       return (tns.tpe().halfCylinderNumber(id) == 1)
0600                  ? (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].first
0601                  : (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].second;
0602     case PclHLS::TPEPanelDisk3:
0603       return (tns.tpe().halfCylinderNumber(id) == 1)
0604                  ? (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].first
0605                  : (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].second;
0606     case PclHLS::TPEPanelDiskM1:
0607       return (tns.tpe().halfCylinderNumber(id) == 1)
0608                  ? (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].first
0609                  : (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].second;
0610     case PclHLS::TPEPanelDiskM2:
0611       return (tns.tpe().halfCylinderNumber(id) == 1)
0612                  ? (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].first
0613                  : (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].second;
0614     case PclHLS::TPEPanelDiskM3:
0615       return (tns.tpe().halfCylinderNumber(id) == 1)
0616                  ? (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].first
0617                  : (tns.tpe().bladeNumber(id) * 2 - (tns.tpe().panelNumber(id) % 2)) + indexHelper[HLS].second;
0618     default:
0619       return -200;
0620   }
0621 }
0622 
0623 //=============================================================================
0624 //===   STATIC CONST MEMBER DEFINITION                                      ===
0625 //=============================================================================
0626 constexpr std::array<double, 6> MillePedeFileReader::multiplier_;
0627 
0628 void MillePedeFileReader::fillPSetDescription(edm::ParameterSetDescription& desc) {
0629   desc.add<std::string>("fileDir", std::string());
0630   desc.add<bool>("ignoreInactiveAlignables", true);
0631   desc.add<std::string>("millePedeEndFile", "millepede.end");
0632   desc.add<std::string>("millePedeLogFile", "millepede.log");
0633   desc.add<std::string>("millePedeResFile", "millepede.res");
0634   desc.add<bool>("isHG", false);
0635 }