CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 <memory>
7 #include <fstream>
8 #include <string>
9 
10 /*** core framework functionality ***/
13 
14 
15 
17 
18  //========================== PUBLIC METHODS ==================================
19  public: //====================================================================
20 
21  explicit MillePedeFileReader(const edm::ParameterSet&);
23 
24  void read();
25  bool storeAlignments();
26 
27  std::array<double, 6> const& getXobs() const { return Xobs; }
28  std::array<double, 6> const& getXobsErr() const { return XobsErr; }
29  std::array<double, 6> const& getTXobs() const { return tXobs; }
30  std::array<double, 6> const& getTXobsErr() const { return tXobsErr; }
31 
32  std::array<double, 6> const& getYobs() const { return Yobs; }
33  std::array<double, 6> const& getYobsErr() const { return YobsErr; }
34  std::array<double, 6> const& getTYobs() const { return tYobs; }
35  std::array<double, 6> const& getTYobsErr() const { return tYobsErr; }
36 
37  std::array<double, 6> const& getZobs() const { return Zobs; }
38  std::array<double, 6> const& getZobsErr() const { return ZobsErr; }
39  std::array<double, 6> const& getTZobs() const { return tZobs; }
40  std::array<double, 6> const& getTZobsErr() const { return tZobsErr; }
41 
42  //========================= PRIVATE METHODS ==================================
43  private: //===================================================================
44 
45  void readMillePedeLogFile();
47 
48  //========================== PRIVATE DATA ====================================
49  //============================================================================
50 
51  // file-names
54 
55  // signifiance of movement must be above
56  double sigCut_;
57  // cutoff in micro-meter & micro-rad
58  double Xcut_, tXcut_;
59  double Ycut_, tYcut_;
60  double Zcut_, tZcut_;
61  // maximum movement in micro-meter/rad
63 
64  double Cutoffs[6] = { Xcut_, Ycut_, Zcut_,
65  tXcut_, tYcut_, tZcut_};
66 
67  bool PedeSuccess = false;
68  bool Movements = false;
69  bool Error = false;
70  bool Significant = false;
71  bool updateDB = false;
72  bool HitMax = false;
73  bool HitErrorMax = false;
74 
75  int Nrec = 0;
76 
77 
78 
79  std::array<double, 6> Xobs = {{0.,0.,0.,0.,0.,0.}};
80  std::array<double, 6> XobsErr = {{0.,0.,0.,0.,0.,0.}};
81  std::array<double, 6> tXobs = {{0.,0.,0.,0.,0.,0.}};
82  std::array<double, 6> tXobsErr = {{0.,0.,0.,0.,0.,0.}};
83 
84  std::array<double, 6> Yobs = {{0.,0.,0.,0.,0.,0.}};
85  std::array<double, 6> YobsErr = {{0.,0.,0.,0.,0.,0.}};
86  std::array<double, 6> tYobs = {{0.,0.,0.,0.,0.,0.}};
87  std::array<double, 6> tYobsErr = {{0.,0.,0.,0.,0.,0.}};
88 
89  std::array<double, 6> Zobs = {{0.,0.,0.,0.,0.,0.}};
90  std::array<double, 6> ZobsErr = {{0.,0.,0.,0.,0.,0.}};
91  std::array<double, 6> tZobs = {{0.,0.,0.,0.,0.,0.}};
92  std::array<double, 6> tZobsErr = {{0.,0.,0.,0.,0.,0.}};
93 
94 };
95 
96 
97 
98 //=============================================================================
99 //=== PUBLIC METHOD IMPLEMENTATION ===
100 //=============================================================================
101 
104  millePedeLogFile_(config.getParameter<std::string>("millePedeLogFile")),
105  millePedeResFile_(config.getParameter<std::string>("millePedeResFile")),
106 
107  sigCut_ (config.getParameter<double>("sigCut")),
108  Xcut_ (config.getParameter<double>("Xcut")),
109  tXcut_ (config.getParameter<double>("tXcut")),
110  Ycut_ (config.getParameter<double>("Ycut")),
111  tYcut_ (config.getParameter<double>("tYcut")),
112  Zcut_ (config.getParameter<double>("Zcut")),
113  tZcut_ (config.getParameter<double>("tZcut")),
114  maxMoveCut_ (config.getParameter<double>("maxMoveCut")),
115  maxErrorCut_(config.getParameter<double>("maxErrorCut"))
116 {
117 }
118 
121  readMillePedeLogFile();
122  readMillePedeResultFile();
123 }
124 
127  return updateDB;
128 }
129 
130 
131 
132 //=============================================================================
133 //=== PRIVATE METHOD IMPLEMENTATION ===
134 //=============================================================================
135 
138 {
139  std::ifstream logFile;
140  logFile.open(millePedeLogFile_.c_str());
141 
142  if (logFile.is_open()) {
143  edm::LogInfo("MillePedeFileReader") << "Reading millepede log-file";
145 
146  while (getline(logFile, line)) {
147  std::string Nrec_string = "NREC =";
148 
149  if (line.find(Nrec_string) != std::string::npos) {
150  std::istringstream iss(line);
151  std::string trash;
152  iss >> trash >> trash >> Nrec;
153 
154  if (Nrec < 25000) {
155  PedeSuccess = false;
156  Movements = false;
157  Error = false;
158  Significant = false;
159  updateDB = false;
160  }
161  }
162  }
163 
164  } else {
165  edm::LogError("MillePedeFileReader") << "Could not read millepede log-file.";
166 
167  PedeSuccess = false;
168  Movements = false;
169  Error = false;
170  Significant = false;
171  updateDB = false;
172  Nrec = 0;
173  }
174 }
175 
178 {
179  std::ifstream resFile;
180  resFile.open(millePedeResFile_.c_str());
181 
182  if (resFile.is_open()) {
183  edm::LogInfo("MillePedeFileReader") << "Reading millepede result-file";
184  double Multiplier[6] = {10000.,10000.,10000.,1000000.,1000000.,1000000.};
185 
187  getline(resFile, line); // drop first line
188 
189  while (getline(resFile, line)) {
190  std::istringstream iss(line);
191 
192  std::vector<std::string> tokens;
194  while (iss >> token) {
195  tokens.push_back(token);
196  }
197 
198  if (tokens.size() > 4 /*3*/) {
199  PedeSuccess = true;
200 
201  int alignable = std::stoi(tokens[0]);
202  int alignableIndex = alignable % 10 - 1;
203 
204  double ObsMove = std::stof(tokens[3]) * Multiplier[alignableIndex];
205  double ObsErr = std::stof(tokens[4]) * Multiplier[alignableIndex];
206 
207  int det = -1;
208 
209  if (alignable >= 60 && alignable <= 69) {
210  det = 2; // TPBHalfBarrel (x+)
211  } else if (alignable >= 8780 && alignable <= 8789) {
212  det = 3; // TPBHalfBarrel (x-)
213  } else if (alignable >= 17520 && alignable <= 17529) {
214  det = 4; // TPEHalfCylinder (x+,z+)
215  } else if (alignable >= 22380 && alignable <= 22389) {
216  det = 5; // TPEHalfCylinder (x-,z+)
217  } else if (alignable >= 27260 && alignable <= 27269) {
218  det = 0; // TPEHalfCylinder (x+,z-)
219  } else if (alignable >= 32120 && alignable <= 32129) {
220  det = 1; //TPEHalfCylinder (x-,z-)
221  } else {
222  continue;
223  }
224 
225  if (alignableIndex == 0 && det >= 0 && det <= 5) {
226  Xobs[det] = ObsMove;
227  XobsErr[det] = ObsErr;
228  } else if (alignableIndex == 1 && det >= 0 && det <= 5) {
229  Yobs[det] = ObsMove;
230  YobsErr[det] = ObsErr;
231  } else if (alignableIndex == 2 && det >= 0 && det <= 5) {
232  Zobs[det] = ObsMove;
233  ZobsErr[det] = ObsErr;
234  } else if (alignableIndex == 3 && det >= 0 && det <= 5) {
235  tXobs[det] = ObsMove;
236  tXobsErr[det] = ObsErr;
237  } else if (alignableIndex == 4 && det >= 0 && det <= 5) {
238  tYobs[det] = ObsMove;
239  tYobsErr[det] = ObsErr;
240  } else if (alignableIndex == 5 && det >= 0 && det <= 5) {
241  tZobs[det] = ObsMove;
242  tZobsErr[det] = ObsErr;
243  }
244 
245  if (abs(ObsMove) > maxMoveCut_) {
246  Movements = false;
247  Error = false;
248  Significant = false;
249  updateDB = false;
250  HitMax = false;
251  continue;
252 
253  } else if (abs(ObsMove) > Cutoffs[alignableIndex]) {
254  Movements = true;
255 
256  if (abs(ObsErr) > maxErrorCut_) {
257  Error = false;
258  Significant = false;
259  updateDB = false;
260  HitErrorMax = true;
261  continue;
262  } else {
263  Error = true;
264  if (abs(ObsMove/ObsErr) > sigCut_) {
265  Significant = true;
266  }
267  }
268  }
269  updateDB = true;
270  }
271  }
272  } else {
273  edm::LogError("MillePedeFileReader") << "Could not read millepede result-file.";
274 
275  PedeSuccess = false;
276  Movements = false;
277  Error = false;
278  Significant = false;
279  updateDB = false;
280  Nrec = 0;
281  }
282 }
283 
284 
285 
286 #endif /* ALIGNMENT_MILLEPEDEALIGNMENTALGORITHM_INTERFACE_MILLEPEDEFILEREADER_H_ */
std::array< double, 6 > const & getTYobs() const
std::array< double, 6 > const & getYobsErr() const
std::array< double, 6 > Yobs
std::array< double, 6 > tZobs
std::array< double, 6 > XobsErr
std::array< double, 6 > Zobs
std::array< double, 6 > tXobsErr
std::array< double, 6 > ZobsErr
std::array< double, 6 > const & getYobs() const
std::array< double, 6 > const & getTXobs() const
std::array< double, 6 > const & getXobsErr() const
std::array< double, 6 > const & getTZobs() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::array< double, 6 > tZobsErr
std::array< double, 6 > tYobs
std::array< double, 6 > Xobs
std::array< double, 6 > tYobsErr
std::array< double, 6 > const & getTXobsErr() const
std::array< double, 6 > const & getZobs() const
std::array< double, 6 > YobsErr
tuple logFile
Definition: heppy_check.py:39
std::array< double, 6 > tXobs
std::array< double, 6 > const & getZobsErr() const
MillePedeFileReader(const edm::ParameterSet &)
std::array< double, 6 > const & getTZobsErr() const
std::array< double, 6 > const & getTYobsErr() const
std::array< double, 6 > const & getXobs() const