CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes
MillePedeFileReader Class Reference

#include <MillePedeFileReader.h>

Public Member Functions

const int binariesAmount () const
 
const mpPCLresults getResults () const
 
const
AlignPCLThresholds::threshold_map 
getThresholdMap () const
 
const std::array< double, 6 > & getTXobs () const
 
const std::array< double, 6 > & getTXobsErr () const
 
const std::array< double, 6 > & getTYobs () const
 
const std::array< double, 6 > & getTYobsErr () const
 
const std::array< double, 6 > & getTZobs () const
 
const std::array< double, 6 > & getTZobsErr () const
 
const std::array< double, 6 > & getXobs () const
 
const std::array< double, 6 > & getXobsErr () const
 
const std::array< double, 6 > & getYobs () const
 
const std::array< double, 6 > & getYobsErr () const
 
const std::array< double, 6 > & getZobs () const
 
const std::array< double, 6 > & getZobsErr () const
 
 MillePedeFileReader (const edm::ParameterSet &, const std::shared_ptr< const PedeLabelerBase > &, const std::shared_ptr< const AlignPCLThresholds > &)
 
void read ()
 
bool storeAlignments ()
 
virtual ~MillePedeFileReader ()=default
 

Private Types

enum  PclHLS : int {
  PclHLS::NotInPCL = -1, PclHLS::TPBHalfBarrelXplus = 2, PclHLS::TPBHalfBarrelXminus = 3, PclHLS::TPEHalfCylinderXplusZplus = 4,
  PclHLS::TPEHalfCylinderXminusZplus = 5, PclHLS::TPEHalfCylinderXplusZminus = 0, PclHLS::TPEHalfCylinderXminusZminus = 1
}
 

Private Member Functions

PclHLS getHLS (const Alignable *)
 
std::string getStringFromHLS (PclHLS HLS)
 
void readMillePedeEndFile ()
 
void readMillePedeLogFile ()
 
void readMillePedeResultFile ()
 

Private Attributes

int binariesAmount_ {0}
 
int exitCode_ {-1}
 
std::string exitMessage_ {""}
 
const std::string millePedeEndFile_
 
const std::string millePedeLogFile_
 
const std::string millePedeResFile_
 
int Nrec_ {0}
 
const std::shared_ptr< const
PedeLabelerBase
pedeLabeler_
 
const std::shared_ptr< const
AlignPCLThresholds
theThresholds_
 
std::array< double, 6 > tXobs_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > tXobsErr_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > tYobs_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > tYobsErr_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > tZobs_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > tZobsErr_ = {{0., 0., 0., 0., 0., 0.}}
 
std::bitset< 4 > updateBits_
 
bool updateDB_ {false}
 
bool vetoUpdateDB_ {false}
 
std::array< double, 6 > Xobs_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > XobsErr_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > Yobs_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > YobsErr_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > Zobs_ = {{0., 0., 0., 0., 0., 0.}}
 
std::array< double, 6 > ZobsErr_ = {{0., 0., 0., 0., 0., 0.}}
 

Static Private Attributes

static constexpr std::array
< double, 6 > 
multiplier_
 

Detailed Description

Definition at line 57 of file MillePedeFileReader.h.

Member Enumeration Documentation

enum MillePedeFileReader::PclHLS : int
strongprivate
Enumerator
NotInPCL 
TPBHalfBarrelXplus 
TPBHalfBarrelXminus 
TPEHalfCylinderXplusZplus 
TPEHalfCylinderXminusZplus 
TPEHalfCylinderXplusZminus 
TPEHalfCylinderXminusZminus 

Definition at line 96 of file MillePedeFileReader.h.

96  : int {
97  NotInPCL = -1,
98  TPBHalfBarrelXplus = 2,
99  TPBHalfBarrelXminus = 3,
100  TPEHalfCylinderXplusZplus = 4,
101  TPEHalfCylinderXminusZplus = 5,
102  TPEHalfCylinderXplusZminus = 0,
103  TPEHalfCylinderXminusZminus = 1
104  };

Constructor & Destructor Documentation

MillePedeFileReader::MillePedeFileReader ( const edm::ParameterSet config,
const std::shared_ptr< const PedeLabelerBase > &  pedeLabeler,
const std::shared_ptr< const AlignPCLThresholds > &  theThresholds 
)
explicit

Definition at line 15 of file MillePedeFileReader.cc.

18  : pedeLabeler_(pedeLabeler),
19  theThresholds_(theThresholds),
20  millePedeEndFile_(config.getParameter<std::string>("millePedeEndFile")),
21  millePedeLogFile_(config.getParameter<std::string>("millePedeLogFile")),
22  millePedeResFile_(config.getParameter<std::string>("millePedeResFile")) {}
const std::shared_ptr< const AlignPCLThresholds > theThresholds_
const std::string millePedeLogFile_
const std::string millePedeResFile_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const std::shared_ptr< const PedeLabelerBase > pedeLabeler_
const std::string millePedeEndFile_
virtual MillePedeFileReader::~MillePedeFileReader ( )
virtualdefault

Member Function Documentation

const int MillePedeFileReader::binariesAmount ( ) const
inline

Definition at line 86 of file MillePedeFileReader.h.

References binariesAmount_.

MillePedeFileReader::PclHLS MillePedeFileReader::getHLS ( const Alignable alignable)
private

Definition at line 235 of file MillePedeFileReader.cc.

References Alignable::alignableObjectId(), Exception, Alignable::id(), NotInPCL, pedeLabeler_, edm::tns(), align::TPBHalfBarrel, TPBHalfBarrelXminus, TPBHalfBarrelXplus, align::TPEHalfCylinder, TPEHalfCylinderXminusZminus, TPEHalfCylinderXminusZplus, TPEHalfCylinderXplusZminus, and TPEHalfCylinderXplusZplus.

Referenced by readMillePedeResultFile().

235  {
236  if (!alignable)
237  return PclHLS::NotInPCL;
238 
239  const auto& tns = pedeLabeler_->alignableTracker()->trackerNameSpace();
240 
241  switch (alignable->alignableObjectId()) {
243  switch (tns.tpb().halfBarrelNumber(alignable->id())) {
244  case 1:
246  case 2:
248  default:
249  throw cms::Exception("LogicError") << "@SUB=MillePedeFileReader::getHLS\n"
250  << "Found a pixel half-barrel number that should not exist: "
251  << tns.tpb().halfBarrelNumber(alignable->id());
252  }
254  switch (tns.tpe().endcapNumber(alignable->id())) {
255  case 1:
256  switch (tns.tpe().halfCylinderNumber(alignable->id())) {
257  case 1:
259  case 2:
261  default:
262  throw cms::Exception("LogicError") << "@SUB=MillePedeFileReader::getHLS\n"
263  << "Found a pixel half-cylinder number that should not exist: "
264  << tns.tpe().halfCylinderNumber(alignable->id());
265  }
266  case 2:
267  switch (tns.tpe().halfCylinderNumber(alignable->id())) {
268  case 1:
270  case 2:
272  default:
273  throw cms::Exception("LogicError") << "@SUB=MillePedeFileReader::getHLS\n"
274  << "Found a pixel half-cylinder number that should not exist: "
275  << tns.tpe().halfCylinderNumber(alignable->id());
276  }
277  default:
278  throw cms::Exception("LogicError")
279  << "@SUB=MillePedeFileReader::getHLS\n"
280  << "Found a pixel endcap number that should not exist: " << tns.tpe().endcapNumber(alignable->id());
281  }
282  default:
283  return PclHLS::NotInPCL;
284  }
285 }
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:180
pathNames_ & tns()), endPathNames_(&tns.getEndPaths()), wantSummary_(tns.wantSummary()
Definition: Schedule.cc:691
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
const std::shared_ptr< const PedeLabelerBase > pedeLabeler_
const mpPCLresults MillePedeFileReader::getResults ( ) const
inline
std::string MillePedeFileReader::getStringFromHLS ( MillePedeFileReader::PclHLS  HLS)
private

Definition at line 287 of file MillePedeFileReader.cc.

References Exception, TPBHalfBarrelXminus, TPBHalfBarrelXplus, TPEHalfCylinderXminusZminus, TPEHalfCylinderXminusZplus, TPEHalfCylinderXplusZminus, and TPEHalfCylinderXplusZplus.

Referenced by readMillePedeResultFile().

287  {
288  switch (HLS) {
290  return "TPBHalfBarrelXminus";
292  return "TPBHalfBarrelXplus";
294  return "TPEHalfCylinderXminusZminus";
296  return "TPEHalfCylinderXplusZminus";
298  return "TPEHalfCylinderXminusZplus";
300  return "TPEHalfCylinderXplusZplus";
301  default:
302  throw cms::Exception("LogicError")
303  << "@SUB=MillePedeFileReader::getStringFromHLS\n"
304  << "Found an alignable structure not possible to map in the default AlignPCLThresholds partitions";
305  }
306 }
const AlignPCLThresholds::threshold_map MillePedeFileReader::getThresholdMap ( ) const
inline

Definition at line 84 of file MillePedeFileReader.h.

References theThresholds_.

84 { return theThresholds_.get()->getThreshold_Map(); }
const std::shared_ptr< const AlignPCLThresholds > theThresholds_
const std::array<double, 6>& MillePedeFileReader::getTXobs ( ) const
inline

Definition at line 71 of file MillePedeFileReader.h.

References tXobs_.

71 { return tXobs_; }
std::array< double, 6 > tXobs_
const std::array<double, 6>& MillePedeFileReader::getTXobsErr ( ) const
inline

Definition at line 72 of file MillePedeFileReader.h.

References tXobsErr_.

72 { return tXobsErr_; }
std::array< double, 6 > tXobsErr_
const std::array<double, 6>& MillePedeFileReader::getTYobs ( ) const
inline

Definition at line 76 of file MillePedeFileReader.h.

References tYobs_.

76 { return tYobs_; }
std::array< double, 6 > tYobs_
const std::array<double, 6>& MillePedeFileReader::getTYobsErr ( ) const
inline

Definition at line 77 of file MillePedeFileReader.h.

References tYobsErr_.

77 { return tYobsErr_; }
std::array< double, 6 > tYobsErr_
const std::array<double, 6>& MillePedeFileReader::getTZobs ( ) const
inline

Definition at line 81 of file MillePedeFileReader.h.

References tZobs_.

81 { return tZobs_; }
std::array< double, 6 > tZobs_
const std::array<double, 6>& MillePedeFileReader::getTZobsErr ( ) const
inline

Definition at line 82 of file MillePedeFileReader.h.

References tZobsErr_.

82 { return tZobsErr_; }
std::array< double, 6 > tZobsErr_
const std::array<double, 6>& MillePedeFileReader::getXobs ( ) const
inline

Definition at line 69 of file MillePedeFileReader.h.

References Xobs_.

69 { return Xobs_; }
std::array< double, 6 > Xobs_
const std::array<double, 6>& MillePedeFileReader::getXobsErr ( ) const
inline

Definition at line 70 of file MillePedeFileReader.h.

References XobsErr_.

70 { return XobsErr_; }
std::array< double, 6 > XobsErr_
const std::array<double, 6>& MillePedeFileReader::getYobs ( ) const
inline

Definition at line 74 of file MillePedeFileReader.h.

References Yobs_.

74 { return Yobs_; }
std::array< double, 6 > Yobs_
const std::array<double, 6>& MillePedeFileReader::getYobsErr ( ) const
inline

Definition at line 75 of file MillePedeFileReader.h.

References YobsErr_.

75 { return YobsErr_; }
std::array< double, 6 > YobsErr_
const std::array<double, 6>& MillePedeFileReader::getZobs ( ) const
inline

Definition at line 79 of file MillePedeFileReader.h.

References Zobs_.

79 { return Zobs_; }
std::array< double, 6 > Zobs_
const std::array<double, 6>& MillePedeFileReader::getZobsErr ( ) const
inline

Definition at line 80 of file MillePedeFileReader.h.

References ZobsErr_.

80 { return ZobsErr_; }
std::array< double, 6 > ZobsErr_
void MillePedeFileReader::read ( )
void MillePedeFileReader::readMillePedeEndFile ( )
private

Definition at line 35 of file MillePedeFileReader.cc.

References exitCode_, exitMessage_, geometryCSVtoXML::line, millePedeEndFile_, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by read().

35  {
36  std::ifstream endFile;
37  endFile.open(millePedeEndFile_.c_str());
38 
39  if (endFile.is_open()) {
40  edm::LogInfo("MillePedeFileReader") << "Reading millepede end-file";
42  getline(endFile, line);
43  std::string trash;
44  if (line.find("-1") != std::string::npos) {
45  getline(endFile, line);
47  std::istringstream iss(line);
48  iss >> exitCode_ >> trash;
49  edm::LogInfo("MillePedeFileReader")
50  << " Pede exit code is: " << exitCode_ << " (" << exitMessage_ << ")" << std::endl;
51  } else {
53  std::istringstream iss(line);
54  iss >> exitCode_ >> trash;
55  edm::LogInfo("MillePedeFileReader")
56  << " Pede exit code is: " << exitCode_ << " (" << exitMessage_ << ")" << std::endl;
57  }
58  } else {
59  edm::LogError("MillePedeFileReader") << "Could not read millepede end-file.";
60  exitMessage_ = "no exit code found";
61  }
62 }
Log< level::Error, false > LogError
Log< level::Info, false > LogInfo
const std::string millePedeEndFile_
void MillePedeFileReader::readMillePedeLogFile ( )
private

Definition at line 64 of file MillePedeFileReader.cc.

References binariesAmount_, geometryCSVtoXML::line, heppy_check::logFile, millePedeLogFile_, Nrec_, AlCaHLTBitMon_QueryRunRegistry::string, theThresholds_, and updateDB_.

Referenced by read().

64  {
65  std::ifstream logFile;
66  logFile.open(millePedeLogFile_.c_str());
67 
68  if (logFile.is_open()) {
69  edm::LogInfo("MillePedeFileReader") << "Reading millepede log-file";
71 
72  while (getline(logFile, line)) {
73  std::string Nrec_string = "NREC =";
74  std::string Binaries_string = "C_binary";
75 
76  if (line.find(Nrec_string) != std::string::npos) {
77  std::istringstream iss(line);
78  std::string trash;
79  iss >> trash >> trash >> Nrec_;
80 
81  if (Nrec_ < theThresholds_->getNrecords()) {
82  edm::LogInfo("MillePedeFileReader")
83  << "Number of records used " << theThresholds_->getNrecords() << std::endl;
84  updateDB_ = false;
85  }
86  }
87 
88  if (line.find(Binaries_string) != std::string::npos) {
89  binariesAmount_ += 1;
90  }
91  }
92  } else {
93  edm::LogError("MillePedeFileReader") << "Could not read millepede log-file.";
94 
95  updateDB_ = false;
96  Nrec_ = 0;
97  }
98 }
Log< level::Error, false > LogError
const std::shared_ptr< const AlignPCLThresholds > theThresholds_
const std::string millePedeLogFile_
Log< level::Info, false > LogInfo
tuple logFile
Definition: heppy_check.py:40
void MillePedeFileReader::readMillePedeResultFile ( )
private

Definition at line 100 of file MillePedeFileReader.cc.

References funct::abs(), getHLS(), getStringFromHLS(), geometryCSVtoXML::line, millePedeResFile_, multiplier_, NotInPCL, Nrec_, pedeLabeler_, AlCaHLTBitMon_QueryRunRegistry::string, AlignPCLThresholds::theta_X, AlignPCLThresholds::theta_Y, AlignPCLThresholds::theta_Z, theThresholds_, unpackBuffers-CaloStage2::token, tXobs_, tXobsErr_, tYobs_, tYobsErr_, tZobs_, tZobsErr_, updateBits_, updateDB_, vetoUpdateDB_, AlignPCLThresholds::X, Xobs_, XobsErr_, AlignPCLThresholds::Y, Yobs_, YobsErr_, AlignPCLThresholds::Z, Zobs_, and ZobsErr_.

Referenced by read().

100  {
101  // cutoffs by coordinate and by alignable
102  std::map<std::string, std::array<float, 6> > cutoffs_;
103  std::map<std::string, std::array<float, 6> > significances_;
104  std::map<std::string, std::array<float, 6> > thresholds_;
105  std::map<std::string, std::array<float, 6> > errors_;
106 
107  std::vector<std::string> alignables_ = theThresholds_->getAlignableList();
108  for (auto& ali : alignables_) {
109  cutoffs_[ali] = theThresholds_->getCut(ali);
110  significances_[ali] = theThresholds_->getSigCut(ali);
111  thresholds_[ali] = theThresholds_->getMaxMoveCut(ali);
112  errors_[ali] = theThresholds_->getMaxErrorCut(ali);
113  }
114 
115  updateDB_ = false;
116  vetoUpdateDB_ = false;
117  std::ifstream resFile;
118  resFile.open(millePedeResFile_.c_str());
119 
120  if (resFile.is_open()) {
121  edm::LogInfo("MillePedeFileReader") << "Reading millepede result-file";
122 
124  getline(resFile, line); // drop first line
125 
126  while (getline(resFile, line)) {
127  std::istringstream iss(line);
128 
129  std::vector<std::string> tokens;
131  while (iss >> token) {
132  tokens.push_back(token);
133  }
134 
135  if (tokens.size() > 4 /*3*/) {
136  auto alignableLabel = std::stoul(tokens[0]);
137  auto alignableIndex = alignableLabel % 10 - 1;
138  const auto alignable = pedeLabeler_->alignableFromLabel(alignableLabel);
139 
140  double ObsMove = std::stof(tokens[3]) * multiplier_[alignableIndex];
141  double ObsErr = std::stof(tokens[4]) * multiplier_[alignableIndex];
142 
143  auto det = getHLS(alignable);
144  int detIndex = static_cast<int>(det);
145  auto coord = static_cast<AlignPCLThresholds::coordType>(alignableIndex);
146  std::string detLabel = getStringFromHLS(det);
147 
148  if (det != PclHLS::NotInPCL) {
149  switch (coord) {
151  Xobs_[detIndex] = ObsMove;
152  XobsErr_[detIndex] = ObsErr;
153  break;
155  Yobs_[detIndex] = ObsMove;
156  YobsErr_[detIndex] = ObsErr;
157  break;
159  Zobs_[detIndex] = ObsMove;
160  ZobsErr_[detIndex] = ObsErr;
161  break;
163  tXobs_[detIndex] = ObsMove;
164  tXobsErr_[detIndex] = ObsErr;
165  break;
167  tYobs_[detIndex] = ObsMove;
168  tYobsErr_[detIndex] = ObsErr;
169  break;
171  tZobs_[detIndex] = ObsMove;
172  tZobsErr_[detIndex] = ObsErr;
173  break;
174  default:
175  edm::LogError("MillePedeFileReader") << "Currently not able to handle DOF " << coord << std::endl;
176  break;
177  }
178  } else {
179  continue;
180  }
181 
182  edm::LogVerbatim("MillePedeFileReader")
183  << " alignableLabel: " << alignableLabel << " with alignableIndex " << alignableIndex << " detIndex"
184  << detIndex << "\n"
185  << " i.e. detLabel: " << detLabel << " (" << coord << ")\n"
186  << " has movement: " << ObsMove << " +/- " << ObsErr << "\n"
187  << " cutoff (cutoffs_[" << detLabel << "][" << coord << "]): " << cutoffs_[detLabel][alignableIndex] << "\n"
188  << " significance (significances_[" << detLabel << "][" << coord
189  << "]): " << significances_[detLabel][alignableIndex] << "\n"
190  << " error thresolds (errors_[" << detLabel << "][" << coord << "]): " << errors_[detLabel][alignableIndex]
191  << "\n"
192  << " max movement (thresholds_[" << detLabel << "][" << coord
193  << "]): " << thresholds_[detLabel][alignableIndex] << "\n"
194  << "=============" << std::endl;
195 
196  if (std::abs(ObsMove) > thresholds_[detLabel][alignableIndex]) {
197  edm::LogWarning("MillePedeFileReader") << "Aborting payload creation."
198  << " Exceeding maximum thresholds for movement: " << std::abs(ObsMove)
199  << " for" << detLabel << "(" << coord << ")";
200  updateBits_.set(0);
201  vetoUpdateDB_ = true;
202  continue;
203 
204  } else if (std::abs(ObsMove) > cutoffs_[detLabel][alignableIndex]) {
205  updateBits_.set(1);
206 
207  if (std::abs(ObsErr) > errors_[detLabel][alignableIndex]) {
208  edm::LogWarning("MillePedeFileReader") << "Aborting payload creation."
209  << " Exceeding maximum thresholds for error: " << std::abs(ObsErr)
210  << " for" << detLabel << "(" << coord << ")";
211  updateBits_.set(2);
212  vetoUpdateDB_ = true;
213  continue;
214  } else {
215  if (std::abs(ObsMove / ObsErr) < significances_[detLabel][alignableIndex]) {
216  updateBits_.set(3);
217  continue;
218  }
219  }
220  updateDB_ = true;
221  edm::LogInfo("MillePedeFileReader")
222  << "This correction: " << ObsMove << "+/-" << ObsErr << " for " << detLabel << "(" << coord
223  << ") will trigger a new Tracker Alignment payload!";
224  }
225  }
226  }
227  } else {
228  edm::LogError("MillePedeFileReader") << "Could not read millepede result-file.";
229 
230  updateDB_ = false;
231  Nrec_ = 0;
232  }
233 }
Log< level::Info, true > LogVerbatim
std::array< double, 6 > tYobs_
std::array< double, 6 > tXobs_
std::bitset< 4 > updateBits_
std::array< double, 6 > YobsErr_
PclHLS getHLS(const Alignable *)
Log< level::Error, false > LogError
static constexpr std::array< double, 6 > multiplier_
const std::shared_ptr< const AlignPCLThresholds > theThresholds_
std::array< double, 6 > tZobs_
std::array< double, 6 > tXobsErr_
std::array< double, 6 > tYobsErr_
std::array< double, 6 > Yobs_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::array< double, 6 > tZobsErr_
std::string getStringFromHLS(PclHLS HLS)
Log< level::Info, false > LogInfo
const std::string millePedeResFile_
std::array< double, 6 > Zobs_
std::array< double, 6 > Xobs_
const std::shared_ptr< const PedeLabelerBase > pedeLabeler_
std::array< double, 6 > ZobsErr_
Log< level::Warning, false > LogWarning
std::array< double, 6 > XobsErr_
bool MillePedeFileReader::storeAlignments ( )

Member Data Documentation

int MillePedeFileReader::binariesAmount_ {0}
private

Definition at line 148 of file MillePedeFileReader.h.

Referenced by binariesAmount(), and readMillePedeLogFile().

int MillePedeFileReader::exitCode_ {-1}
private

Definition at line 151 of file MillePedeFileReader.h.

Referenced by getResults(), and readMillePedeEndFile().

std::string MillePedeFileReader::exitMessage_ {""}
private

Definition at line 152 of file MillePedeFileReader.h.

Referenced by getResults(), and readMillePedeEndFile().

const std::string MillePedeFileReader::millePedeEndFile_
private

Definition at line 125 of file MillePedeFileReader.h.

Referenced by readMillePedeEndFile().

const std::string MillePedeFileReader::millePedeLogFile_
private

Definition at line 126 of file MillePedeFileReader.h.

Referenced by readMillePedeLogFile().

const std::string MillePedeFileReader::millePedeResFile_
private

Definition at line 127 of file MillePedeFileReader.h.

Referenced by readMillePedeResultFile().

constexpr std::array< double, 6 > MillePedeFileReader::multiplier_
staticprivate
Initial value:
= {{10000.,
10000.,
10000.,
1000000.,
1000000.,
1000000.}}

Definition at line 130 of file MillePedeFileReader.h.

Referenced by readMillePedeResultFile().

int MillePedeFileReader::Nrec_ {0}
private
const std::shared_ptr<const PedeLabelerBase> MillePedeFileReader::pedeLabeler_
private

Definition at line 119 of file MillePedeFileReader.h.

Referenced by getHLS(), and readMillePedeResultFile().

const std::shared_ptr<const AlignPCLThresholds> MillePedeFileReader::theThresholds_
private
std::array<double, 6> MillePedeFileReader::tXobs_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 156 of file MillePedeFileReader.h.

Referenced by getTXobs(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::tXobsErr_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 157 of file MillePedeFileReader.h.

Referenced by getTXobsErr(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::tYobs_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 161 of file MillePedeFileReader.h.

Referenced by getTYobs(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::tYobsErr_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 162 of file MillePedeFileReader.h.

Referenced by getTYobsErr(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::tZobs_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 166 of file MillePedeFileReader.h.

Referenced by getTZobs(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::tZobsErr_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 167 of file MillePedeFileReader.h.

Referenced by getTZobsErr(), and readMillePedeResultFile().

std::bitset<4> MillePedeFileReader::updateBits_
private

Definition at line 145 of file MillePedeFileReader.h.

Referenced by getResults(), and readMillePedeResultFile().

bool MillePedeFileReader::updateDB_ {false}
private
bool MillePedeFileReader::vetoUpdateDB_ {false}
private

Definition at line 138 of file MillePedeFileReader.h.

Referenced by getResults(), readMillePedeResultFile(), and storeAlignments().

std::array<double, 6> MillePedeFileReader::Xobs_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 154 of file MillePedeFileReader.h.

Referenced by getXobs(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::XobsErr_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 155 of file MillePedeFileReader.h.

Referenced by getXobsErr(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::Yobs_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 159 of file MillePedeFileReader.h.

Referenced by getYobs(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::YobsErr_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 160 of file MillePedeFileReader.h.

Referenced by getYobsErr(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::Zobs_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 164 of file MillePedeFileReader.h.

Referenced by getZobs(), and readMillePedeResultFile().

std::array<double, 6> MillePedeFileReader::ZobsErr_ = {{0., 0., 0., 0., 0., 0.}}
private

Definition at line 165 of file MillePedeFileReader.h.

Referenced by getZobsErr(), and readMillePedeResultFile().