CMS 3D CMS Logo

ME0GeometryValidate.cc
Go to the documentation of this file.
1 /*
2 //\class ME0GeometryValidate
3 
4  Description: ME0 GeometryValidate for DD4hep
5 
6  //
7 // Author: Sergio Lo Meo (sergio.lo.meo@cern.ch) following what Ianna Osborne made for DTs (DD4hep migration)
8 // Created: 29 Apr 2020
9 */
10 
16 
22 
24 
27 
28 #include <TFile.h>
29 #include <TH1.h>
30 
31 #include <limits>
32 #include <string>
33 #include <type_traits>
34 #include <algorithm>
35 
36 using namespace std;
37 
38 template <class T>
39 typename enable_if<!numeric_limits<T>::is_integer, bool>::type almost_equal(T x, T y, int ulp) {
40  // the machine epsilon has to be scaled to the magnitude of the values used
41  // and multiplied by the desired precision in ULPs (units in the last place)
42  return abs(x - y) <= numeric_limits<T>::epsilon() * abs(x + y) * ulp
43  // unless the result is subnormal
44  || abs(x - y) < numeric_limits<T>::min();
45 }
46 
47 using namespace edm;
48 
49 class ME0GeometryValidate : public one::EDAnalyzer<> {
50 public:
51  explicit ME0GeometryValidate(const ParameterSet&);
52  ~ME0GeometryValidate() override {}
53 
54 private:
55  void beginJob() override;
56  void analyze(const edm::Event&, const edm::EventSetup&) override;
57  void endJob() override;
58 
59  void validateME0ChamberGeometry();
60  void validateME0EtaPartitionGeometry();
61  void validateME0EtaPartitionGeometry2();
62 
63  void compareTransform(const GlobalPoint&, const TGeoMatrix*);
64  void compareShape(const GeomDet*, const float*);
65 
66  float getDistance(const GlobalPoint&, const GlobalPoint&);
67  float getDiff(const float, const float);
68 
69  void makeHistograms(const char*);
70  void makeHistograms2(const char*);
71  void makeHistogram(const string&, vector<float>&);
72 
73  void clearData() {
74  globalDistances_.clear();
75  topWidths_.clear();
76  bottomWidths_.clear();
77  lengths_.clear();
78  thicknesses_.clear();
79  }
80 
81  void clearData2() {
82  nstrips_.clear();
83  pitch_.clear();
84  stripslen_.clear();
85  }
86 
90  TFile* outFile_;
91  vector<float> globalDistances_;
92  vector<float> topWidths_;
93  vector<float> bottomWidths_;
94  vector<float> lengths_;
95  vector<float> thicknesses_;
96  vector<float> nstrips_;
97  vector<float> pitch_;
98  vector<float> stripslen_;
99  string infileName_;
100  string outfileName_;
102 };
103 
105  : tokGeom_{esConsumes<ME0Geometry, MuonGeometryRecord>(edm::ESInputTag{})},
106  infileName_(iConfig.getUntrackedParameter<string>("infileName", "cmsRecoGeom-2026.root")),
107  outfileName_(iConfig.getUntrackedParameter<string>("outfileName", "validateME0Geometry.root")),
108  tolerance_(iConfig.getUntrackedParameter<int>("tolerance", 6)) {
109  fwGeometry_.loadMap(infileName_.c_str());
110  outFile_ = TFile::Open(outfileName_.c_str(), "RECREATE");
111 }
112 
114  me0Geometry_ = &eventSetup.getData(tokGeom_);
115  LogTrace("ME0Geometry") << "Validating ME0 chamber geometry";
119 }
120 
122  clearData();
123 
124  for (auto const& it : me0Geometry_->chambers()) {
125  ME0DetId chId = it->id();
126  GlobalPoint gp = it->surface().toGlobal(LocalPoint(0.0, 0.0, 0.0));
127  const TGeoMatrix* matrix = fwGeometry_.getMatrix(chId.rawId());
128 
129  if (!matrix) {
130  LogTrace("ME0Geometry") << "Failed to get matrix of ME0 chamber with detid: " << chId.rawId();
131  continue;
132  }
134 
135  auto const& shape = fwGeometry_.getShapePars(chId.rawId());
136 
137  if (!shape) {
138  LogTrace("ME0Geometry") << "Failed to get shape of ME0 chamber with detid: " << chId.rawId();
139  continue;
140  }
142  }
143  makeHistograms("ME0 Chamber");
144 }
145 
147  clearData();
148 
149  for (auto const& it : me0Geometry_->etaPartitions()) {
150  ME0DetId chId = it->id();
151  GlobalPoint gp = it->surface().toGlobal(LocalPoint(0.0, 0.0, 0.0));
152  const TGeoMatrix* matrix = fwGeometry_.getMatrix(chId.rawId());
153 
154  if (!matrix) {
155  LogTrace("ME0Geometry") << "Failed to get matrix of ME0 eta partition 2 with detid: " << chId.rawId();
156  continue;
157  }
159 
160  auto const& shape = fwGeometry_.getShapePars(chId.rawId());
161 
162  if (!shape) {
163  LogTrace("ME0Geometry") << "Failed to get shape of ME0 eta partition 2 with detid: " << chId.rawId();
164  continue;
165  }
167  }
168  makeHistograms("ME0 Eta Partition");
169 }
170 
172  clearData2();
173 
174  for (auto const& it : me0Geometry_->etaPartitions()) {
175  ME0DetId chId = it->id();
176  const int n_strips = it->nstrips();
177  const float n_pitch = it->pitch();
178  const StripTopology& topo = it->specificTopology();
179  const float stripLen = topo.stripLength();
180  const float* parameters = fwGeometry_.getParameters(chId.rawId());
181  nstrips_.push_back(std::abs(n_strips - parameters[0]));
182 
183  if (n_strips) {
184  for (int istrips = 1; istrips <= n_strips; istrips++) {
185  if (std::abs(n_pitch - parameters[2]) < 1.0e-05)
186  pitch_.push_back(0.0f);
187  else
188  pitch_.push_back(std::abs(n_pitch - parameters[2]));
189  if (std::abs(stripLen - parameters[1]) < 1.0e-05)
190  pitch_.push_back(0.0f);
191  else
192  stripslen_.push_back(std::abs(stripLen - parameters[1]));
193  }
194  } else {
195  LogWarning("ME0Geometry") << "ATTENTION! nStrips == 0";
196  }
197  }
198  makeHistograms2("ME0 Eta Partition");
199 }
200 
202  double local[3] = {0.0, 0.0, 0.0};
203  double global[3];
204 
205  matrix->LocalToMaster(local, global);
206 
207  float distance = getDistance(GlobalPoint(global[0], global[1], global[2]), gp);
208  if (abs(distance) < 1.0e-7)
209  distance = 0.0; // set a tollerance for the distance inside Histos
210  globalDistances_.push_back(distance);
211 }
212 
213 void ME0GeometryValidate::compareShape(const GeomDet* det, const float* shape) {
214  float shapeTopWidth;
215  float shapeBottomWidth;
216  float shapeLength;
217  float shapeThickness;
218 
219  if (shape[0] == 1) {
220  shapeTopWidth = shape[2];
221  shapeBottomWidth = shape[1];
222  shapeLength = shape[4];
223  shapeThickness = shape[3];
224  } else if (shape[0] == 2) {
225  shapeTopWidth = shape[1];
226  shapeBottomWidth = shape[1];
227  shapeLength = shape[2];
228  shapeThickness = shape[3];
229  } else {
230  LogTrace("ME0Geometry") << "Failed to get box or trapezoid from shape";
231 
232  return;
233  }
234 
235  float topWidth, bottomWidth;
236  float length, thickness;
237 
238  const Bounds* bounds = &(det->surface().bounds());
239  if (const TrapezoidalPlaneBounds* tpbs = dynamic_cast<const TrapezoidalPlaneBounds*>(bounds)) {
240  array<const float, 4> const& ps = tpbs->parameters();
241 
242  assert(ps.size() == 4);
243 
244  bottomWidth = ps[0];
245  topWidth = ps[1];
246  thickness = ps[2];
247  length = ps[3];
248 
249  } else if ((dynamic_cast<const RectangularPlaneBounds*>(bounds))) {
250  length = det->surface().bounds().length() * 0.5;
251  topWidth = det->surface().bounds().width() * 0.5;
252  bottomWidth = topWidth;
253  thickness = det->surface().bounds().thickness() * 0.5;
254 
255  } else {
256  LogTrace("ME0Geometry") << "Failed to get bounds";
257 
258  return;
259  }
260  topWidths_.push_back(std::abs(shapeTopWidth - topWidth));
261  bottomWidths_.push_back(std::abs(shapeBottomWidth - bottomWidth));
262  lengths_.push_back(std::abs(shapeLength - length));
263  thicknesses_.push_back(std::abs(shapeThickness - thickness));
264 }
265 
267  return sqrt((p1.x() - p2.x()) * (p1.x() - p2.x()) + (p1.y() - p2.y()) * (p1.y() - p2.y()) +
268  (p1.z() - p2.z()) * (p1.z() - p2.z()));
269 }
270 
271 float ME0GeometryValidate::getDiff(const float val1, const float val2) {
272  if (almost_equal(val1, val2, tolerance_))
273  return 0.0f;
274  else
275  return (val1 - val2);
276 }
277 
279  outFile_->cd();
280 
281  string d(detector);
282 
283  string gdn = d + ": distance between points in global coordinates";
285 
286  string twn = d + ": absolute difference between top widths (along X)";
288 
289  string bwn = d + ": absolute difference between bottom widths (along X)";
291 
292  string ln = d + ": absolute difference between lengths (along Y)";
293  makeHistogram(ln, lengths_);
294 
295  string tn = d + ": absolute difference between thicknesses (along Z)";
297 }
298 
300  outFile_->cd();
301 
302  string d(detector);
303 
304  string ns = d + ": absolute difference between nStrips";
305  makeHistogram(ns, nstrips_);
306 
307  string pi = d + ": absolute difference between Strips Pitch";
309 
310  string pl = d + ": absolute difference between Strips Length";
312 }
313 
314 void ME0GeometryValidate::makeHistogram(const string& name, vector<float>& data) {
315  if (data.empty())
316  return;
317 
318  const auto [minE, maxE] = minmax_element(begin(data), end(data));
319 
320  TH1D hist(name.c_str(), name.c_str(), 100, *minE * (1 + 0.10), *maxE * (1 + 0.10));
321 
322  for (auto const& it : data)
323  hist.Fill(it);
324 
325  hist.GetXaxis()->SetTitle("[cm]");
326  hist.Write();
327 }
328 
330 
332  LogTrace("ME0Geometry") << "Done.";
333  LogTrace("ME0Geometry") << "Results written to " << outfileName_;
334  outFile_->Close();
335 }
336 
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
virtual float length() const =0
float getDiff(const float, const float)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
enable_if<!numeric_limits< T >::is_integer, bool >::type almost_equal(T x, T y, int ulp)
assert(be >=bs)
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
#define LogTrace(id)
void beginJob()
Definition: Breakpoints.cc:14
const Double_t pi
vector< float > thicknesses_
void analyze(const edm::Event &, const edm::EventSetup &) override
virtual float thickness() const =0
vector< float > nstrips_
virtual float stripLength() const =0
vector< float > globalDistances_
const std::vector< const ME0Chamber * > & chambers() const
Return a vector of all ME0 chambers.
Definition: ME0Geometry.cc:29
T sqrt(T t)
Definition: SSEVec.h:19
const ME0Geometry * me0Geometry_
const edm::ESGetToken< ME0Geometry, MuonGeometryRecord > tokGeom_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float getDistance(const GlobalPoint &, const GlobalPoint &)
vector< float > topWidths_
d
Definition: ztail.py:151
vector< float > lengths_
void compareShape(const GeomDet *, const float *)
void compareTransform(const GlobalPoint &, const TGeoMatrix *)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
void makeHistograms2(const char *)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const float * getShapePars(unsigned int id) const
Definition: FWGeometry.cc:461
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:450
const TGeoMatrix * getMatrix(unsigned int id) const
Definition: FWGeometry.cc:225
float x
void makeHistograms(const char *)
vector< float > stripslen_
vector< float > bottomWidths_
void makeHistogram(const string &, vector< float > &)
Definition: Bounds.h:18
Log< level::Warning, false > LogWarning
long double T
virtual float width() const =0
ME0GeometryValidate(const ParameterSet &)
const std::vector< ME0EtaPartition const * > & etaPartitions() const
Return a vector of all ME0 eta partitions.
Definition: ME0Geometry.cc:33
Definition: event.py:1
const Bounds & bounds() const
Definition: Surface.h:87