CMS 3D CMS Logo

GeometryInterface.h
Go to the documentation of this file.
1 #ifndef SiPixel_GeometryInterface_h
2 #define SiPixel_GeometryInterface_h
3 // -*- C++ -*-
4 //
5 // Package: SiPixelPhase1Common
6 // Class: GeometryInterface
7 //
8 // The histogram manager uses this class to gather information about a sample.
9 // All geometry dependence goes here.
10 //
11 // Original Author: Marcel Schneider
12 //
13 
18 
19 #include <functional>
20 #include <map>
21 #include <string>
22 #include <array>
23 
25  public:
26  // an ID is produced by interning a string name.
27  typedef int ID;
28  // A column could have multiple IDs if it is a or-form.
29  // Not used atm, makes many things much easier.
30  typedef ID Column;
31  typedef double Value;
32  static const Value UNDEFINED;
33 
34  // Essentially a map backed by a vector (for the small counts here
35  // this should always be faster). Most ops turned out to be not needed.
36  typedef std::vector<std::pair<Column, Value>> Values;
37 
39 
40  bool loaded() { return is_loaded; };
41 
42  // The hard work happens here.
43  void load(edm::EventSetup const& iSetup);
44 
46  // in this order the struct should fit 2 64bit words and is cheap to copy.
49  int16_t col = 0;
50  int16_t row = 0;
51  };
52 
53  // This has to be fast, _should_ not malloc.
54  void extractColumns(std::vector<Column> const& names,
55  InterestingQuantities const& iq, Values& out) {
56  out.clear();
57  for (Column const& col : names) {
58  auto val = extract(col, iq);
59  out.push_back(val);
60  }
61  }
62 
63  // the pair return value is historical; it is only really needed with or-columns.
64  // But it is cleaner to carry it around.
65  std::pair<Column, Value> extract(Column const& col, InterestingQuantities const& iq) {
66  assert(col != 0 || !"Extracting invalid column.");
67  ID id = col;
68  assert(ID(extractors.size()) > id || !"extractors vector too small!");
69  auto& ex = extractors[id];
70  if (!ex) { // we have never heard about this. This is a typo for sure.
71  edm::LogError("GeometryInterface")
72  << "Undefined column used: " << unintern(id)
73  << ". Check your spelling.\n";
74  } else {
75  auto val = ex(iq);
76  if (val != UNDEFINED) {
77  return std::make_pair(Column{id}, val);
78  }
79  }
80  return std::make_pair(col, UNDEFINED);
81  }
82 
83  Value extract(ID id, DetId did, edm::Event* ev = nullptr, int16_t col = 0,
84  int16_t row = 0) {
85  InterestingQuantities iq = {ev, did, col, row};
86  return extractors[id](iq);
87  }
88 
89  // TODO: for Phase0 (and maybe also Phase2) this should include the 4 corners
90  // of each ROC (or the *correct* corners of the respective modules).
91  std::vector<InterestingQuantities> const& allModules() {
92  return all_modules;
93  }
94 
95  Value maxValue(ID id) { return max_value[id]; };
96  Value minValue(ID id) { return min_value[id]; };
97  Value binWidth(ID id) { return bin_width[id]; };
98 
99  // turn string into an ID, adding it if needed.
100  ID intern(std::string const& id) {
101  auto it = ids.find(id);
102  if (it == ids.end()) {
103  ids[id] = ++max_id;
104  extractors.resize(max_id + 1);
105  }
106  return ids[id];
107  };
108 
109  // turn an ID back into a string. Only for pretty output (including histo
110  // labels), so it can be slow (though intern() does not have to be fast
111  // either).
113  for (auto& e : ids)
114  if (e.second == id) return e.first;
115  return "INVALID";
116  }
117 
119  return unintern(col);
120  }
121 
122  std::string formatValue(Column, Value);
123 
124  private:
125  void loadFromTopology(edm::EventSetup const& iSetup,
126  const edm::ParameterSet& iConfig);
127  void loadFromSiPixelCoordinates(edm::EventSetup const& iSetup,
128  const edm::ParameterSet& iConfig);
129  void loadTimebased(edm::EventSetup const& iSetup,
130  const edm::ParameterSet& iConfig);
131  void loadModuleLevel(edm::EventSetup const& iSetup,
132  const edm::ParameterSet& iConfig);
133  void loadFEDCabling(edm::EventSetup const& iSetup,
134  const edm::ParameterSet& iConfig);
135 
137 
138  bool is_loaded = false;
139 
140  // This holds closures that compute the column values in step1.
141  // can be a Vector since ids are dense.
142  std::vector<std::function<Value(InterestingQuantities const& iq)>> extractors;
143  // quantity range if it is known. Can be UNDEFINED, in this case booking will
144  // determine the range. Map for ease of use.
145  std::map<ID, Value> max_value;
146  std::map<ID, Value> min_value;
147  std::map<ID, Value> bin_width;
148 
149  // cache of pre-formatted values. Can be pre-populated while loading
150  // (used for Pixel*Name)
151  std::map<std::pair<Column, Value>, std::string> format_value;
152 
153  void addExtractor(ID id,
155  Value min = UNDEFINED, Value max = UNDEFINED, Value binwidth = 1) {
156  max_value[id] = max;
157  min_value[id] = min;
158  bin_width[id] = binwidth;
159  extractors[id] = func;
160  }
161 
162  std::vector<InterestingQuantities> all_modules;
163 
164  // interning table. Maps string IDs to a dense set of integer IDs
165  std::map<std::string, ID> ids{std::make_pair(std::string("INVALID"), ID(0))};
166  ID max_id = 0;
167 };
168 
169 #endif
std::map< ID, Value > max_value
std::map< std::pair< Column, Value >, std::string > format_value
std::map< std::string, ID > ids
static const HistoName names[]
GeometryInterface(const edm::ParameterSet &conf)
std::string pretty(Column col)
void addExtractor(ID id, std::function< Value(InterestingQuantities const &iq)> func, Value min=UNDEFINED, Value max=UNDEFINED, Value binwidth=1)
std::string formatValue(Column, Value)
Value binWidth(ID id)
bool ev
std::string unintern(ID id)
void loadFromTopology(edm::EventSetup const &iSetup, const edm::ParameterSet &iConfig)
Value maxValue(ID id)
void load(edm::EventSetup const &iSetup)
void extractColumns(std::vector< Column > const &names, InterestingQuantities const &iq, Values &out)
const edm::ParameterSet iConfig
std::vector< InterestingQuantities > all_modules
std::map< ID, Value > bin_width
T min(T a, T b)
Definition: MathUtil.h:58
static const Value UNDEFINED
std::vector< std::pair< Column, Value > > Values
Value extract(ID id, DetId did, edm::Event *ev=0, int16_t col=0, int16_t row=0)
Definition: DetId.h:18
void loadFEDCabling(edm::EventSetup const &iSetup, const edm::ParameterSet &iConfig)
Value minValue(ID id)
std::pair< Column, Value > extract(Column const &col, InterestingQuantities const &iq)
void loadModuleLevel(edm::EventSetup const &iSetup, const edm::ParameterSet &iConfig)
std::vector< std::function< Value(InterestingQuantities const &iq)> > extractors
std::vector< InterestingQuantities > const & allModules()
std::map< ID, Value > min_value
ID intern(std::string const &id)
void loadTimebased(edm::EventSetup const &iSetup, const edm::ParameterSet &iConfig)
void loadFromSiPixelCoordinates(edm::EventSetup const &iSetup, const edm::ParameterSet &iConfig)