CMS 3D CMS Logo

printStripTrackerMap.cc
Go to the documentation of this file.
2 #include <cstdlib>
3 #include <fstream>
4 #include <iostream>
5 #include <numeric> // std::accumulate
6 #include <sstream>
7 #include <string>
8 #include <vector>
9 
10 #include "TCanvas.h"
11 #include "TStyle.h"
12 
13 int main(int argc, char* argv[]) {
15  std::vector<std::pair<uint32_t, float>> detidValues;
16 
17  // Parse command line arguments
18  for (int i = 1; i < argc; ++i) {
19  if (std::string(argv[i]) == "--input-file" && i + 1 < argc) {
20  gStyle->SetPalette(kRainbow);
21  gStyle->SetNumberContours(256);
22  inputFile = argv[++i];
23  } else {
24  gStyle->SetPalette(1);
25  // Treat as DetId list if no --input-file is provided
26  uint32_t detid = std::stoul(argv[i]);
27  detidValues.emplace_back(detid, 1.0); // Default value is 1.0
28  }
29  }
30 
31  // If --input-file is provided, read from file
32  if (!inputFile.empty()) {
33  std::ifstream file(inputFile);
34  if (!file) {
35  std::cerr << "Error: Unable to open input file " << inputFile << std::endl;
36  return 1;
37  }
38 
40  while (std::getline(file, line)) {
41  std::istringstream iss(line);
42  uint32_t detid;
43  float value = 1.0; // Default value
44 
45  iss >> detid;
46  if (iss >> value) { // If a second column exists, read it as value
47  detidValues.emplace_back(detid, value);
48  } else {
49  detidValues.emplace_back(detid, 1.0);
50  }
51  }
52  }
53 
54  // Create the map and fill it
55  SiStripTkMaps theMap("COLZ0 AL");
56  theMap.bookMap("Strip Tracker Map of Marked modules", "input values");
57 
58  for (const auto& [detid, value] : detidValues) {
59  theMap.fill(detid, value);
60  }
61 
62  // Check if all values are the same using a lambda function
63  bool allSame = std::all_of(detidValues.begin(), detidValues.end(), [&](const std::pair<uint32_t, float>& p) {
64  return p.second == detidValues[0].second;
65  });
66 
67  TCanvas c = TCanvas("c", "c");
68  theMap.drawMap(c, "");
69 
70  // adjust the z-axis scale
71  if (allSame)
72  theMap.setZAxisRange(0., detidValues[0].second);
73 
74  c.SaveAs("SiStripsTkMaps.png");
75 
76  std::cout << "Filled tracker map with " << detidValues.size() << " detids." << std::endl;
77 
78  return 0;
79 }
int main(int argc, char *argv[])
void bookMap(const std::string mapTitle, const std::string zAxisTitle)
U second(std::pair< T, U > const &p)
Definition: value.py:1
void fill(long rawid, double val)
void setZAxisRange(double xmin, double xmax)
Definition: SiStripTkMaps.h:69
void drawMap(TCanvas &canvas, std::string option="")