CMS 3D CMS Logo

MillePedeFileReader.h
Go to the documentation of this file.
1 #ifndef ALIGNMENT_MILLEPEDEALIGNMENTALGORITHM_INTERFACE_MILLEPEDEFILEREADER_H_
2 #define ALIGNMENT_MILLEPEDEALIGNMENTALGORITHM_INTERFACE_MILLEPEDEFILEREADER_H_
3 
4 /*** system includes ***/
5 #include <array>
6 #include <string>
7 #include <iostream>
8 
9 /*** core framework functionality ***/
12 
13 /*** Alignment ***/
16 
17 struct mpPCLresults {
18 private:
24  std::bitset<4> m_updateBits;
25 
26 public:
27  mpPCLresults(bool isDBUpdated,
28  bool isDBUpdateVetoed,
29  int nRecords,
30  int exitCode,
31  std::string exitMessage,
32  std::bitset<4> updateBits)
33  : m_isDBUpdated(isDBUpdated),
34  m_isDBUpdateVetoed(isDBUpdateVetoed),
35  m_nRecords(nRecords),
36  m_exitCode(exitCode),
37  m_exitMessage(exitMessage),
38  m_updateBits(updateBits) {}
39 
40  const bool getDBUpdated() { return m_isDBUpdated; }
41  const bool getDBVetoed() { return m_isDBUpdateVetoed; }
42  const bool exceedsThresholds() { return m_updateBits.test(0); }
43  const bool exceedsCutoffs() { return m_updateBits.test(1); }
44  const bool exceedsMaxError() { return m_updateBits.test(2); }
45  const bool belowSignificance() { return m_updateBits.test(3); }
46  const int getNRecords() { return m_nRecords; }
47  const int getExitCode() { return m_exitCode; }
49 
50  void print() {
51  edm::LogInfo("MillePedeFileReader") << " is DB updated: " << m_isDBUpdated
52  << " is DB update vetoed: " << m_isDBUpdateVetoed << " nRecords: " << m_nRecords
53  << " exitCode: " << m_exitCode << " (" << m_exitMessage << ")" << std::endl;
54  }
55 };
56 
58  //========================== PUBLIC METHODS ==================================
59 public: //====================================================================
60  explicit MillePedeFileReader(const edm::ParameterSet&,
61  const std::shared_ptr<const PedeLabelerBase>&,
62  const std::shared_ptr<const AlignPCLThresholds>&);
63 
64  virtual ~MillePedeFileReader() = default;
65 
66  void read();
67  bool storeAlignments();
68 
69  const std::array<double, 6>& getXobs() const { return Xobs_; }
70  const std::array<double, 6>& getXobsErr() const { return XobsErr_; }
71  const std::array<double, 6>& getTXobs() const { return tXobs_; }
72  const std::array<double, 6>& getTXobsErr() const { return tXobsErr_; }
73 
74  const std::array<double, 6>& getYobs() const { return Yobs_; }
75  const std::array<double, 6>& getYobsErr() const { return YobsErr_; }
76  const std::array<double, 6>& getTYobs() const { return tYobs_; }
77  const std::array<double, 6>& getTYobsErr() const { return tYobsErr_; }
78 
79  const std::array<double, 6>& getZobs() const { return Zobs_; }
80  const std::array<double, 6>& getZobsErr() const { return ZobsErr_; }
81  const std::array<double, 6>& getTZobs() const { return tZobs_; }
82  const std::array<double, 6>& getTZobsErr() const { return tZobsErr_; }
83 
84  const AlignPCLThresholds::threshold_map getThresholdMap() const { return theThresholds_.get()->getThreshold_Map(); }
85 
86  const int binariesAmount() const { return binariesAmount_; }
87 
88  const mpPCLresults getResults() const {
90  }
91 
92 private:
93  //========================= PRIVATE ENUMS ====================================
94  //============================================================================
95 
96  enum class PclHLS : int {
97  NotInPCL = -1,
104  };
105 
106  //========================= PRIVATE METHODS ==================================
107  //============================================================================
108 
109  void readMillePedeEndFile();
110  void readMillePedeLogFile();
112  PclHLS getHLS(const Alignable*);
114 
115  //========================== PRIVATE DATA ====================================
116  //============================================================================
117 
118  // pede labeler plugin
119  const std::shared_ptr<const PedeLabelerBase> pedeLabeler_;
120 
121  // thresholds from DB
122  const std::shared_ptr<const AlignPCLThresholds> theThresholds_;
123 
124  // file-names
128 
129  // conversion factors: cm to um & rad to urad
130  static constexpr std::array<double, 6> multiplier_ = {{10000., // X
131  10000., // Y
132  10000., // Z
133  1000000., // tX
134  1000000., // tY
135  1000000.}}; // tZ
136 
137  bool updateDB_{false};
138  bool vetoUpdateDB_{false};
139 
140  // stores in a compact format the 4 decisions:
141  // 1st bit: exceeds maximum thresholds
142  // 2nd bit: exceeds cutoffs (significant movement)
143  // 3rd bit: exceeds maximum errors
144  // 4th bit: is below the significance
145  std::bitset<4> updateBits_;
146 
147  // pede binaries available
149 
150  int Nrec_{0};
151  int exitCode_{-1};
153 
154  std::array<double, 6> Xobs_ = {{0., 0., 0., 0., 0., 0.}};
155  std::array<double, 6> XobsErr_ = {{0., 0., 0., 0., 0., 0.}};
156  std::array<double, 6> tXobs_ = {{0., 0., 0., 0., 0., 0.}};
157  std::array<double, 6> tXobsErr_ = {{0., 0., 0., 0., 0., 0.}};
158 
159  std::array<double, 6> Yobs_ = {{0., 0., 0., 0., 0., 0.}};
160  std::array<double, 6> YobsErr_ = {{0., 0., 0., 0., 0., 0.}};
161  std::array<double, 6> tYobs_ = {{0., 0., 0., 0., 0., 0.}};
162  std::array<double, 6> tYobsErr_ = {{0., 0., 0., 0., 0., 0.}};
163 
164  std::array<double, 6> Zobs_ = {{0., 0., 0., 0., 0., 0.}};
165  std::array<double, 6> ZobsErr_ = {{0., 0., 0., 0., 0., 0.}};
166  std::array<double, 6> tZobs_ = {{0., 0., 0., 0., 0., 0.}};
167  std::array<double, 6> tZobsErr_ = {{0., 0., 0., 0., 0., 0.}};
168 };
169 
170 const std::array<std::string, 8> coord_str = {{"X", "Y", "Z", "theta_X", "theta_Y", "theta_Z", "extra_DOF", "none"}};
171 inline std::ostream& operator<<(std::ostream& os, const AlignPCLThresholds::coordType& c) {
173  return os << "unrecongnized coordinate";
174  return os << coord_str[c];
175 }
176 
177 #endif /* ALIGNMENT_MILLEPEDEALIGNMENTALGORITHM_INTERFACE_MILLEPEDEFILEREADER_H_ */
MillePedeFileReader::readMillePedeLogFile
void readMillePedeLogFile()
Definition: MillePedeFileReader.cc:64
MillePedeFileReader::getTZobsErr
const std::array< double, 6 > & getTZobsErr() const
Definition: MillePedeFileReader.h:82
MillePedeFileReader::PclHLS
PclHLS
Definition: MillePedeFileReader.h:96
MillePedeFileReader::storeAlignments
bool storeAlignments()
Definition: MillePedeFileReader.cc:30
mpPCLresults
Definition: MillePedeFileReader.h:17
MillePedeFileReader::getTZobs
const std::array< double, 6 > & getTZobs() const
Definition: MillePedeFileReader.h:81
MillePedeFileReader::getXobs
const std::array< double, 6 > & getXobs() const
Definition: MillePedeFileReader.h:69
MessageLogger.h
MillePedeFileReader::binariesAmount_
int binariesAmount_
Definition: MillePedeFileReader.h:148
MillePedeFileReader::multiplier_
static constexpr std::array< double, 6 > multiplier_
Definition: MillePedeFileReader.h:130
MillePedeFileReader::read
void read()
Definition: MillePedeFileReader.cc:24
MillePedeFileReader::updateBits_
std::bitset< 4 > updateBits_
Definition: MillePedeFileReader.h:145
Alignable
Definition: Alignable.h:27
MillePedeFileReader::PclHLS::TPEHalfCylinderXplusZplus
operator<<
std::ostream & operator<<(std::ostream &os, const AlignPCLThresholds::coordType &c)
Definition: MillePedeFileReader.h:171
MillePedeFileReader::Nrec_
int Nrec_
Definition: MillePedeFileReader.h:150
mpPCLresults::getExitCode
const int getExitCode()
Definition: MillePedeFileReader.h:47
mpPCLresults::m_isDBUpdateVetoed
bool m_isDBUpdateVetoed
Definition: MillePedeFileReader.h:20
AlignPCLThresholds::X
Definition: AlignPCLThresholds.h:14
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
MillePedeFileReader::tZobs_
std::array< double, 6 > tZobs_
Definition: MillePedeFileReader.h:166
PedeLabelerBase.h
MillePedeFileReader::getTYobsErr
const std::array< double, 6 > & getTYobsErr() const
Definition: MillePedeFileReader.h:77
MillePedeFileReader::exitCode_
int exitCode_
Definition: MillePedeFileReader.h:151
MillePedeFileReader::getYobs
const std::array< double, 6 > & getYobs() const
Definition: MillePedeFileReader.h:74
MillePedeFileReader::~MillePedeFileReader
virtual ~MillePedeFileReader()=default
MillePedeFileReader::ZobsErr_
std::array< double, 6 > ZobsErr_
Definition: MillePedeFileReader.h:165
MillePedeFileReader::getTYobs
const std::array< double, 6 > & getTYobs() const
Definition: MillePedeFileReader.h:76
mpPCLresults::m_isDBUpdated
bool m_isDBUpdated
Definition: MillePedeFileReader.h:19
MillePedeFileReader::pedeLabeler_
const std::shared_ptr< const PedeLabelerBase > pedeLabeler_
Definition: MillePedeFileReader.h:119
AlignPCLThresholds::endOfTypes
Definition: AlignPCLThresholds.h:14
mpPCLresults::getDBUpdated
const bool getDBUpdated()
Definition: MillePedeFileReader.h:40
mpPCLresults::exceedsMaxError
const bool exceedsMaxError()
Definition: MillePedeFileReader.h:44
MillePedeFileReader::getTXobs
const std::array< double, 6 > & getTXobs() const
Definition: MillePedeFileReader.h:71
MillePedeFileReader::PclHLS::NotInPCL
MillePedeFileReader::millePedeLogFile_
const std::string millePedeLogFile_
Definition: MillePedeFileReader.h:126
MillePedeFileReader::getZobs
const std::array< double, 6 > & getZobs() const
Definition: MillePedeFileReader.h:79
MillePedeFileReader::getXobsErr
const std::array< double, 6 > & getXobsErr() const
Definition: MillePedeFileReader.h:70
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
MillePedeFileReader::getYobsErr
const std::array< double, 6 > & getYobsErr() const
Definition: MillePedeFileReader.h:75
MillePedeFileReader::tYobsErr_
std::array< double, 6 > tYobsErr_
Definition: MillePedeFileReader.h:162
MillePedeFileReader::PclHLS::TPEHalfCylinderXminusZplus
edm::ParameterSet
Definition: ParameterSet.h:47
AlignPCLThresholds.h
MillePedeFileReader::binariesAmount
const int binariesAmount() const
Definition: MillePedeFileReader.h:86
MillePedeFileReader::tYobs_
std::array< double, 6 > tYobs_
Definition: MillePedeFileReader.h:161
mpPCLresults::exceedsCutoffs
const bool exceedsCutoffs()
Definition: MillePedeFileReader.h:43
MillePedeFileReader::readMillePedeResultFile
void readMillePedeResultFile()
Definition: MillePedeFileReader.cc:100
MillePedeFileReader::tZobsErr_
std::array< double, 6 > tZobsErr_
Definition: MillePedeFileReader.h:167
MillePedeFileReader::PclHLS::TPEHalfCylinderXminusZminus
MillePedeFileReader::getZobsErr
const std::array< double, 6 > & getZobsErr() const
Definition: MillePedeFileReader.h:80
mpPCLresults::m_updateBits
std::bitset< 4 > m_updateBits
Definition: MillePedeFileReader.h:24
MillePedeFileReader
Definition: MillePedeFileReader.h:57
MillePedeFileReader::PclHLS::TPBHalfBarrelXminus
mpPCLresults::getExitMessage
const std::string getExitMessage()
Definition: MillePedeFileReader.h:48
createfilelist.int
int
Definition: createfilelist.py:10
MillePedeFileReader::millePedeEndFile_
const std::string millePedeEndFile_
Definition: MillePedeFileReader.h:125
MillePedeFileReader::PclHLS::TPEHalfCylinderXplusZminus
MillePedeFileReader::updateDB_
bool updateDB_
Definition: MillePedeFileReader.h:137
MillePedeFileReader::millePedeResFile_
const std::string millePedeResFile_
Definition: MillePedeFileReader.h:127
MillePedeFileReader::YobsErr_
std::array< double, 6 > YobsErr_
Definition: MillePedeFileReader.h:160
MillePedeFileReader::tXobsErr_
std::array< double, 6 > tXobsErr_
Definition: MillePedeFileReader.h:157
MillePedeFileReader::getThresholdMap
const AlignPCLThresholds::threshold_map getThresholdMap() const
Definition: MillePedeFileReader.h:84
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
AlignPCLThresholds::threshold_map
std::map< std::string, AlignPCLThreshold > threshold_map
Definition: AlignPCLThresholds.h:13
mpPCLresults::belowSignificance
const bool belowSignificance()
Definition: MillePedeFileReader.h:45
MillePedeFileReader::MillePedeFileReader
MillePedeFileReader(const edm::ParameterSet &, const std::shared_ptr< const PedeLabelerBase > &, const std::shared_ptr< const AlignPCLThresholds > &)
Definition: MillePedeFileReader.cc:15
mpPCLresults::m_exitMessage
std::string m_exitMessage
Definition: MillePedeFileReader.h:23
MillePedeFileReader::PclHLS::TPBHalfBarrelXplus
MillePedeFileReader::getResults
const mpPCLresults getResults() const
Definition: MillePedeFileReader.h:88
coord_str
const std::array< std::string, 8 > coord_str
Definition: MillePedeFileReader.h:170
MillePedeFileReader::tXobs_
std::array< double, 6 > tXobs_
Definition: MillePedeFileReader.h:156
MillePedeFileReader::getStringFromHLS
std::string getStringFromHLS(PclHLS HLS)
Definition: MillePedeFileReader.cc:287
mpPCLresults::exceedsThresholds
const bool exceedsThresholds()
Definition: MillePedeFileReader.h:42
mpPCLresults::getDBVetoed
const bool getDBVetoed()
Definition: MillePedeFileReader.h:41
MillePedeFileReader::readMillePedeEndFile
void readMillePedeEndFile()
Definition: MillePedeFileReader.cc:35
AlignPCLThresholds::coordType
coordType
Definition: AlignPCLThresholds.h:14
MillePedeFileReader::Xobs_
std::array< double, 6 > Xobs_
Definition: MillePedeFileReader.h:154
MillePedeFileReader::Zobs_
std::array< double, 6 > Zobs_
Definition: MillePedeFileReader.h:164
MillePedeFileReader::getTXobsErr
const std::array< double, 6 > & getTXobsErr() const
Definition: MillePedeFileReader.h:72
MillePedeFileReader::theThresholds_
const std::shared_ptr< const AlignPCLThresholds > theThresholds_
Definition: MillePedeFileReader.h:122
MillePedeFileReader::XobsErr_
std::array< double, 6 > XobsErr_
Definition: MillePedeFileReader.h:155
ParameterSet.h
mpPCLresults::getNRecords
const int getNRecords()
Definition: MillePedeFileReader.h:46
mpPCLresults::mpPCLresults
mpPCLresults(bool isDBUpdated, bool isDBUpdateVetoed, int nRecords, int exitCode, std::string exitMessage, std::bitset< 4 > updateBits)
Definition: MillePedeFileReader.h:27
mpPCLresults::m_nRecords
int m_nRecords
Definition: MillePedeFileReader.h:21
mpPCLresults::m_exitCode
int m_exitCode
Definition: MillePedeFileReader.h:22
mpPCLresults::print
void print()
Definition: MillePedeFileReader.h:50
MillePedeFileReader::Yobs_
std::array< double, 6 > Yobs_
Definition: MillePedeFileReader.h:159
MillePedeFileReader::exitMessage_
std::string exitMessage_
Definition: MillePedeFileReader.h:152
MillePedeFileReader::vetoUpdateDB_
bool vetoUpdateDB_
Definition: MillePedeFileReader.h:138
MillePedeFileReader::getHLS
PclHLS getHLS(const Alignable *)
Definition: MillePedeFileReader.cc:235