CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCalWaferValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Validation/HGCalValidation
4 // Class: HGCalWaferValidation
5 //
20 //
21 // Original Author: Imran Yusuff
22 // Created: Thu, 27 May 2021 19:47:08 GMT
23 //
24 // Additional Author(s): Yulun Miao
25 /*
26  Changelog:
27 
28  Wed, 17 Nov 2021 21:04:10 UTC by Yulun Miao:
29 
30  * Use the l or h preceding the thickness to determine the type of flat file
31  * Unified partial wafer information to HGCalTypes.h to allow cross-compare
32 
33  Tue, 14 Dmr 2021 by Imran Yusuff:
34 
35  * Further tidying the code
36  * Now uses the first line of the flat file to determine flat file type (single number means new format)
37  * Separate out shape and rotation matching operation into functions
38  * Fixed validation for D86 geometry (especially in rotation for type-3 layers)
39 */
40 
41 // system include files
42 #include <memory>
43 
44 // user include files
47 
50 
55 
59 
60 #include <fstream>
61 #include <regex>
62 
63 //
64 // class declaration
65 //
66 
67 // If the analyzer does not use TFileService, please remove
68 // the template argument to the base class so the class inherits
69 // from edm::one::EDAnalyzer<>
70 // This will improve performance in multithreaded jobs.
71 
73 public:
74  explicit HGCalWaferValidation(const edm::ParameterSet&);
75  ~HGCalWaferValidation() override;
76 
77  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
78 
79 private:
80  // This is MessageLogger logging category
81  const std::string logcat = "HGCalWaferValidation";
82 
83  // wafer coordinate is (layer, u, v), used as index for map
84  using WaferCoord = std::tuple<int, int, int>;
85 
86  std::string strWaferCoord(const WaferCoord& coord);
87 
88  // mapping of wafer shape codes: DD / old format -> new format (LD/HD)
89  using WaferShapeMap = std::map<std::string, int>;
90 
91  const WaferShapeMap waferShapeMapDD = {{"F", HGCalTypes::WaferPartialType::WaferFull},
92  {"a", HGCalTypes::WaferPartialType::WaferHalf},
93  {"am", HGCalTypes::WaferPartialType::WaferHalf2},
94  {"b", HGCalTypes::WaferPartialType::WaferFive},
95  {"bm", HGCalTypes::WaferPartialType::WaferFive2},
96  {"c", HGCalTypes::WaferPartialType::WaferThree},
97  {"d", HGCalTypes::WaferPartialType::WaferSemi},
98  {"dm", HGCalTypes::WaferPartialType::WaferSemi2},
99  {"g", HGCalTypes::WaferPartialType::WaferChopTwo},
100  {"gm", HGCalTypes::WaferPartialType::WaferChopTwoM}};
101 
102  const WaferShapeMap waferShapeMapLD = {{"0", HGCalTypes::WaferPartialType::WaferFull},
103  {"1", HGCalTypes::WaferPartialType::WaferHalf},
104  {"2", HGCalTypes::WaferPartialType::WaferHalf},
105  {"3", HGCalTypes::WaferPartialType::WaferSemi},
106  {"4", HGCalTypes::WaferPartialType::WaferSemi},
107  {"5", HGCalTypes::WaferPartialType::WaferFive},
108  {"6", HGCalTypes::WaferPartialType::WaferThree}};
109 
110  const WaferShapeMap waferShapeMapHD = {{"0", HGCalTypes::WaferPartialType::WaferFull},
111  {"1", HGCalTypes::WaferPartialType::WaferHalf2},
112  {"2", HGCalTypes::WaferPartialType::WaferChopTwoM},
113  {"3", HGCalTypes::WaferPartialType::WaferSemi2},
114  {"4", HGCalTypes::WaferPartialType::WaferSemi2},
115  {"5", HGCalTypes::WaferPartialType::WaferFive2}};
116 
117  bool DDFindHGCal(DDCompactView::GraphWalker& walker, std::string targetName);
120  bool isThicknessMatched(const int geoThickClass, const int fileThickness);
121  bool isRotationMatched(
122  const bool isNewFile, const int layer, const int fileShapeCode, const int geoRotCode, const int fileRotCode);
123 
124  void beginJob() override;
125  void analyze(const edm::Event&, const edm::EventSetup&) override;
126  void endJob() override;
127 
128  // ----------member data ---------------------------
129  // module parameters
131 
132  // information from newer flat file header
133  unsigned int layerCount_;
134  std::vector<int> layerTypes_;
135 
136  // struct to hold wafer information from DD in map
137  struct WaferInfo {
139  double x;
140  double y;
142  int rotCode;
143  std::string waferName; // TEMPORARY
144  };
145 
146  // EDM token to access DD
148 
149  // map holding all wafer properties from DD
150  std::map<WaferCoord, struct WaferInfo> waferData_;
151 
152  // boolean map to keep track of unaccounted DD wafers (not in flat file)
153  std::map<WaferCoord, bool> waferValidated_;
154 };
155 
156 //
157 // constants, enums and typedefs
158 //
159 
160 //
161 // static data member definitions
162 //
163 
164 //
165 // constructors and destructor
166 //
168  : geometryFileName_(iConfig.getParameter<edm::FileInPath>("GeometryFileName")) {
169  viewToken_ = esConsumes<DDCompactView, IdealGeometryRecord>();
170  //now do what ever initialization is needed
171 }
172 
174  // do anything here that needs to be done at desctruction time
175  // (e.g. close files, deallocate resources etc.)
176  //
177  // please remove this method altogether if it would be left empty
178 }
179 
180 //
181 // member functions
182 //
183 
184 // convert WaferCoord tuple to string representation (i.e. for printing)
186  std::stringstream ss;
187  ss << "(" << std::get<0>(coord) << "," << std::get<1>(coord) << "," << std::get<2>(coord) << ")";
188  return ss.str();
189 }
190 
191 // ----- find HGCal entry among the DD -----
193  if (walker.current().first.name().name() == targetName) {
194  // target found
195  return true;
196  }
197  if (walker.firstChild()) {
198  do {
199  if (DDFindHGCal(walker, targetName))
200  // target inside child
201  return true;
202  } while (walker.nextSibling());
203  walker.parent();
204  }
205  return false;
206 }
207 
208 // ----- find the next wafer, then process the wafer layer -----
210  if (walker.current().first.name().fullname().rfind("hgcalwafer:", 0) == 0) {
211  // first wafer found. Now process the entire layer of wafers.
212  ProcessWaferLayer(walker);
213  return;
214  }
215  if (walker.firstChild()) {
216  do {
217  DDFindWafers(walker);
218  } while (walker.nextSibling());
219  walker.parent();
220  }
221 }
222 
223 // ----- process the layer of wafers -----
225  static int waferLayer = 0; // layer numbers in DD are assumed to be sequential from 1
226  waferLayer++;
227  edm::LogVerbatim(logcat) << "ProcessWaferLayer: Processing layer " << waferLayer;
228  do {
229  if (walker.current().first.name().fullname().rfind("hgcalwafer:", 0) == 0) {
230  auto wafer = walker.current();
231  const std::string waferName(walker.current().first.name().fullname());
232  //edm::LogVerbatim(logcat) << " " << waferName; // DEBUG: in case wafer name info is needed
233  const int copyNo = wafer.second->copyno();
234  // extract DD layer properties
235  const int waferType = HGCalTypes::getUnpackedType(copyNo);
236  const int waferU = HGCalTypes::getUnpackedU(copyNo);
237  const int waferV = HGCalTypes::getUnpackedV(copyNo);
238  const WaferCoord waferCoord(waferLayer, waferU, waferV); // map index
239  // build struct of DD wafer properties
240  struct WaferInfo waferInfo;
241  waferInfo.waferName = waferName; //TEMPORARY
242  waferInfo.thickClass = waferType;
243  waferInfo.x = wafer.second->translation().x();
244  waferInfo.y = wafer.second->translation().y();
245  const std::string waferNameData =
246  std::regex_replace(waferName,
247  std::regex("(HGCal[EH]E)(Wafer[01])(Fine|Coarse[12])([a-z]*)([0-9]*)"),
248  "$1 $2-$3 $4 $5",
249  std::regex_constants::format_no_copy);
250  std::stringstream ss(waferNameData);
251  std::string EEorHE;
252  std::string typeStr;
253  std::string shapeStr;
254  std::string rotStr;
255  ss >> EEorHE >> typeStr >> shapeStr >> rotStr;
256  // assume rotational symmetry of full-sized wafers
257  if (shapeStr.empty())
258  shapeStr = "F";
259  if (rotStr.empty())
260  rotStr = "0";
261  const int rotCode(std::stoi(rotStr));
262  //edm::LogVerbatim(logcat) << "rotStr " << rotStr << " rotCode " << rotCode;
263 
264  // convert shape code to wafer types defined in HGCalTypes.h
265  waferInfo.shapeCode = waferShapeMapDD.at(shapeStr);
266 
267  waferInfo.rotCode = rotCode;
268  // populate the map
269  waferData_[waferCoord] = waferInfo;
270  waferValidated_[waferCoord] = false;
271  }
272  } while (walker.nextSibling());
273 }
274 
275 // ------------ check if wafer thickness in geo matches in file (true = is a match) ------------
276 bool HGCalWaferValidation::isThicknessMatched(const int geoThickClass, const int fileThickness) {
277  if (geoThickClass == 0 && fileThickness == 120)
278  return true;
279  if (geoThickClass == 1 && fileThickness == 200)
280  return true;
281  if (geoThickClass == 2 && fileThickness == 300)
282  return true;
283  return false;
284 }
285 
286 // ------------ check if wafer rotation in geo matches in file (true = is a match) ------------
288  const bool isNewFile, const int layer, const int fileShapeCode, const int geoRotCode, const int fileRotCode) {
289  if (fileShapeCode != HGCalTypes::WaferPartialType::WaferFull && geoRotCode == fileRotCode)
290  return true;
291  if (fileShapeCode == HGCalTypes::WaferPartialType::WaferFull) {
292  if (isNewFile && layerTypes_[layer - 1] == 3) { // this array is index-0 based
293  if ((geoRotCode + 1) % 2 == fileRotCode % 2)
294  return true;
295  } else {
296  if (geoRotCode % 2 == fileRotCode % 2)
297  return true;
298  }
299  }
300  return false;
301 }
302 
303 // ------------ method called for each event ------------
305  using namespace edm;
306 
307  // Get the CMS DD
308  auto viewH = iSetup.getHandle(viewToken_);
309 
310  if (!viewH.isValid()) {
311  edm::LogPrint(logcat) << "Error obtaining geometry handle!";
312  return;
313  }
314 
315  edm::LogVerbatim(logcat) << "Root is : " << viewH->root();
316  edm::LogVerbatim(logcat) << std::endl;
317 
318  // find HGCalEE
319  auto eeWalker = viewH->walker();
320  const bool eeFound = DDFindHGCal(eeWalker, "HGCalEE");
321  if (eeFound) {
322  edm::LogVerbatim(logcat) << "HGCalEE found!";
323  edm::LogVerbatim(logcat) << "name = " << eeWalker.current().first.name().name();
324  edm::LogVerbatim(logcat) << "fullname = " << eeWalker.current().first.name().fullname();
325  } else {
326  edm::LogPrint(logcat) << "HGCalEE not found!";
327  }
328  edm::LogVerbatim(logcat) << std::endl;
329 
330  // find HGCalHEsil
331  auto hesilWalker = viewH->walker();
332  const bool hesilFound = DDFindHGCal(hesilWalker, "HGCalHEsil");
333  if (hesilFound) {
334  edm::LogVerbatim(logcat) << "HGCalHEsil found!";
335  edm::LogVerbatim(logcat) << "name = " << hesilWalker.current().first.name().name();
336  edm::LogVerbatim(logcat) << "fullname = " << hesilWalker.current().first.name().fullname();
337  } else {
338  edm::LogPrint(logcat) << "HGCalHEsil not found!";
339  }
340  edm::LogVerbatim(logcat) << std::endl;
341 
342  // find HGCalHEmix
343  auto hemixWalker = viewH->walker();
344  const bool hemixFound = DDFindHGCal(hemixWalker, "HGCalHEmix");
345  if (hemixFound) {
346  edm::LogVerbatim(logcat) << "HGCalHEmix found!";
347  edm::LogVerbatim(logcat) << "name = " << hemixWalker.current().first.name().name();
348  edm::LogVerbatim(logcat) << "fullname = " << hemixWalker.current().first.name().fullname();
349  } else {
350  edm::LogPrint(logcat) << "HGCalHEmix not found!";
351  }
352  edm::LogVerbatim(logcat) << std::endl;
353 
354  // give up if no HGCal found at all
355  if (!(eeFound || hesilFound || hemixFound)) {
356  edm::LogPrint(logcat) << "Nothing found. Giving up.";
357  return;
358  }
359 
360  // Now walk the HGCalEE walker to find the first wafer on each layer and process them
361  edm::LogVerbatim(logcat) << "Calling DDFindWafers(eeWalker);";
362  DDFindWafers(eeWalker);
363 
364  // Walk the HGCalHEsilwalker to find the first wafer on each layer and process them
365  edm::LogVerbatim(logcat) << "Calling DDFindWafers(hesilWalker);";
366  DDFindWafers(hesilWalker);
367 
368  // Walk the HGCalHEmix walker to find the first wafer on each layer and process them
369  edm::LogVerbatim(logcat) << "Calling DDFindWafers(hemixWalker);";
370  DDFindWafers(hemixWalker);
371 
372  // Confirm all the DD wafers have been read
373  edm::LogVerbatim(logcat) << "Number of wafers read from DD: " << waferData_.size();
374 
375  // Now open the geometry text file
377  edm::LogVerbatim(logcat) << "Opening geometry text file: " << fileName;
378  std::ifstream geoTxtFile(fileName);
379 
380  if (!geoTxtFile) {
381  edm::LogPrint(logcat) << "Cannot open geometry text file.";
382  return;
383  }
384 
385  // total processed counter
386  int nTotalProcessed = 0;
387 
388  // geometry error counters
389  int nMissing = 0;
390  int nThicknessError = 0;
391  int nPosXError = 0;
392  int nPosYError = 0;
393  int nShapeError = 0;
394  int nRotError = 0;
395  int nUnaccounted = 0;
396 
398 
399  // find out if this file is an old file or a new file
400  std::getline(geoTxtFile, buf);
401  std::stringstream ss(buf);
402  std::vector<std::string> first_tokens;
403  while (ss >> buf)
404  if (!buf.empty())
405  first_tokens.push_back(buf);
406 
407  const bool isNewFile(first_tokens.size() == 1);
408 
409  if (isNewFile) {
410  edm::LogVerbatim(logcat) << "Text file is of newer version.";
411  layerCount_ = std::stoi(buf);
412  std::getline(geoTxtFile, buf);
413  std::stringstream layerTypesSS(buf);
414  while (layerTypesSS >> buf)
415  if (!buf.empty())
416  layerTypes_.push_back(std::stoi(buf));
417  if (layerTypes_.size() != layerCount_)
418  edm::LogWarning(logcat) << "Number of layer types does not tally with layer count.";
419 
420  // TEMP: make sure reading is correct
421  edm::LogVerbatim(logcat) << "layerCount = " << layerCount_;
422  for (unsigned i = 0; i < layerTypes_.size(); i++) {
423  edm::LogVerbatim(logcat) << " layerType " << i + 1 << " = " << layerTypes_[i]; // 1-based
424  }
425  } else {
426  layerCount_ = 0;
427  // rewind back
428  geoTxtFile.clear();
429  geoTxtFile.seekg(0);
430  }
431 
432  // process each line on the text file
433  while (std::getline(geoTxtFile, buf)) {
434  std::stringstream ss(buf);
435  std::vector<std::string> tokens;
436  while (ss >> buf)
437  if (!buf.empty())
438  tokens.push_back(buf);
439  if (tokens.size() != 8)
440  continue;
441 
442  nTotalProcessed++;
443 
444  // extract wafer info from a textfile line
445  const int waferLayer(std::stoi(tokens[0]));
446  const std::string waferShapeStr(tokens[1]);
447  const std::string waferDensityStr(isNewFile ? tokens[2].substr(0, 1) : "");
448  const int waferThickness(isNewFile ? std::stoi(tokens[2].substr(1)) : std::stoi(tokens[2]));
449  const double waferX(std::stod(tokens[3]));
450  const double waferY(std::stod(tokens[4]));
451  const int waferRotCode(std::stoi(tokens[5]));
452  const int waferU(std::stoi(tokens[6]));
453  const int waferV(std::stoi(tokens[7]));
454  const int waferShapeCode(isNewFile ? (waferDensityStr == "l" ? waferShapeMapLD.at(waferShapeStr)
455  : waferDensityStr == "h" ? waferShapeMapHD.at(waferShapeStr)
456  : HGCalTypes::WaferPartialType::WaferOut)
457  : waferShapeMapDD.at(waferShapeStr));
458 
459  // map index for crosschecking with DD
460  const WaferCoord waferCoord(waferLayer, waferU, waferV);
461 
462  // now check for (and report) wafer data disagreements
463 
464  if (waferData_.find(waferCoord) == waferData_.end()) {
465  nMissing++;
466  edm::LogVerbatim(logcat) << "MISSING: " << strWaferCoord(waferCoord);
467  continue;
468  }
469 
470  const struct WaferInfo waferInfo = waferData_[waferCoord];
471  waferValidated_[waferCoord] = true;
472 
473  if (!isThicknessMatched(waferInfo.thickClass, waferThickness)) {
474  nThicknessError++;
475  edm::LogVerbatim(logcat) << "THICKNESS ERROR: " << strWaferCoord(waferCoord);
476  }
477 
478  // it seems that wafer x-coords relative to their layer plane are mirrored...
479  if (fabs(-waferInfo.x - waferX) > 0.015) { // assuming this much tolerance
480  nPosXError++;
481  edm::LogVerbatim(logcat) << "POSITION x ERROR: " << strWaferCoord(waferCoord);
482  }
483 
484  if (fabs(waferInfo.y - waferY) > 0.015) { // assuming this much tolerance
485  nPosYError++;
486  edm::LogVerbatim(logcat) << "POSITION y ERROR: " << strWaferCoord(waferCoord);
487  }
488 
489  if (waferInfo.shapeCode != waferShapeCode || waferShapeCode == HGCalTypes::WaferPartialType::WaferOut) {
490  nShapeError++;
491  edm::LogVerbatim(logcat) << "SHAPE ERROR: " << strWaferCoord(waferCoord) << " ( " << waferInfo.shapeCode
492  << " != " << waferDensityStr << waferShapeCode << " ) name=" << waferInfo.waferName;
493  }
494 
495  if (!isRotationMatched(isNewFile, std::get<0>(waferCoord), waferShapeCode, waferInfo.rotCode, waferRotCode)) {
496  nRotError++;
497  edm::LogVerbatim(logcat) << "ROTATION ERROR: " << strWaferCoord(waferCoord) << " ( " << waferInfo.rotCode
498  << " != " << waferRotCode << " (" << waferShapeCode
499  << ") ) name=" << waferInfo.waferName;
500  }
501  }
502 
503  geoTxtFile.close();
504 
505  // Find unaccounted DD wafers
506  for (auto const& accounted : waferValidated_) {
507  if (!accounted.second) {
508  nUnaccounted++;
509  edm::LogVerbatim(logcat) << "UNACCOUNTED: " << strWaferCoord(accounted.first);
510  }
511  }
512 
513  // Print out error counts
514  edm::LogVerbatim(logcat) << std::endl;
515  edm::LogVerbatim(logcat) << "*** ERROR COUNTS ***";
516  edm::LogVerbatim(logcat) << "Missing : " << nMissing;
517  edm::LogVerbatim(logcat) << "Thickness error : " << nThicknessError;
518  edm::LogVerbatim(logcat) << "Pos-x error : " << nPosXError;
519  edm::LogVerbatim(logcat) << "Pos-y error : " << nPosYError;
520  edm::LogVerbatim(logcat) << "Shape error : " << nShapeError;
521  edm::LogVerbatim(logcat) << "Rotation error : " << nRotError;
522  edm::LogVerbatim(logcat) << "Unaccounted : " << nUnaccounted;
523  edm::LogVerbatim(logcat) << std::endl;
524  edm::LogVerbatim(logcat) << "Total wafers processed from geotxtfile = " << nTotalProcessed;
525  edm::LogVerbatim(logcat) << std::endl;
526 
527  // Issue a LogPrint (warning) if there is at least one wafer errors
528  if (nMissing > 0 || nThicknessError > 0 || nPosXError > 0 || nPosYError > 0 || nShapeError > 0 || nRotError > 0 ||
529  nUnaccounted > 0) {
530  edm::LogPrint(logcat) << "There are at least one wafer error.";
531  }
532 }
533 
534 // ------------ method called once each job just before starting event loop ------------
536  // please remove this method if not needed
537 }
538 
539 // ------------ method called once each job just after ending the event loop ------------
541  // please remove this method if not needed
542 }
543 
544 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
546  //The following says we do not know what parameters are allowed so do no validation
547  // Please change this to state exactly what you do use, even if it is no parameters
549  desc.add<edm::FileInPath>("GeometryFileName",
550  edm::FileInPath("Validation/HGCalValidation/data/geomnew_corrected_360.txt"));
551  descriptions.add("hgcalWaferValidation", desc);
552 }
553 
554 //define this as a plug-in
std::map< WaferCoord, struct WaferInfo > waferData_
Log< level::Info, true > LogVerbatim
const WaferShapeMap waferShapeMapHD
void ProcessWaferLayer(DDCompactView::GraphWalker &walker)
bool DDFindHGCal(DDCompactView::GraphWalker &walker, std::string targetName)
const WaferShapeMap waferShapeMapDD
static int32_t getUnpackedU(int id)
Definition: HGCalTypes.cc:16
int32_t waferU(const int32_t index)
std::map< std::string, int > WaferShapeMap
int32_t waferLayer(const int32_t index)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const std::string logcat
static int32_t getUnpackedV(int id)
Definition: HGCalTypes.cc:22
std::tuple< int, int, int > WaferCoord
constexpr std::array< uint8_t, layerIndexSize > layer
const WaferShapeMap waferShapeMapLD
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
result_type parent()
Definition: GraphWalker.h:130
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< int > layerTypes_
bool isThicknessMatched(const int geoThickClass, const int fileThickness)
static int32_t getUnpackedType(int id)
Definition: HGCalTypes.cc:14
result_type firstChild()
Definition: GraphWalker.h:108
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isRotationMatched(const bool isNewFile, const int layer, const int fileShapeCode, const int geoRotCode, const int fileRotCode)
Log< level::Warning, true > LogPrint
result_type nextSibling()
Definition: GraphWalker.h:119
edm::ESGetToken< DDCompactView, IdealGeometryRecord > viewToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::map< WaferCoord, bool > waferValidated_
std::string strWaferCoord(const WaferCoord &coord)
int32_t waferV(const int32_t index)
std::string fullPath() const
Definition: FileInPath.cc:161
edm::FileInPath geometryFileName_
void DDFindWafers(DDCompactView::GraphWalker &walker)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
HGCalWaferValidation(const edm::ParameterSet &)
value_type current() const
Definition: GraphWalker.h:86