CMS 3D CMS Logo

DTGeometryValidate.cc
Go to the documentation of this file.
7 
11 
13 
16 
17 #include <TFile.h>
18 #include <TH1.h>
19 
20 #include <cmath>
21 #include <limits>
22 #include <string>
23 #include <type_traits>
24 #include <algorithm>
25 
26 using namespace std;
27 
28 template <class T>
29 typename enable_if<!numeric_limits<T>::is_integer, bool>::type almost_equal(T x, T y, int ulp) {
30  // the machine epsilon has to be scaled to the magnitude of the values used
31  // and multiplied by the desired precision in ULPs (units in the last place)
32  return abs(x - y) <= numeric_limits<T>::epsilon() * abs(x + y) * ulp
33  // unless the result is subnormal
34  || abs(x - y) < numeric_limits<T>::min();
35 }
36 
37 using namespace edm;
38 
39 class DTGeometryValidate : public one::EDAnalyzer<> {
40 public:
41  explicit DTGeometryValidate(const ParameterSet&);
42  ~DTGeometryValidate() override {}
43 
44 private:
45  void beginJob() override;
46  void analyze(const edm::Event&, const edm::EventSetup&) override;
47  void endJob() override;
48 
49  void validateDTChamberGeometry();
50  void validateDTLayerGeometry();
51 
52  void compareTransform(const GlobalPoint&, const TGeoMatrix*);
53  void compareShape(const GeomDet*, const float*);
54 
55  float getDistance(const GlobalPoint&, const GlobalPoint&);
56  float getDiff(const float, const float);
57 
58  void makeHistograms(const char*);
59  void makeHistogram(const string&, vector<float>&);
60 
61  void clearData() {
62  globalDistances_.clear();
63  topWidths_.clear();
64  bottomWidths_.clear();
65  lengths_.clear();
66  thicknesses_.clear();
67  }
68 
72  TFile* outFile_;
73  vector<float> globalDistances_;
74  vector<float> topWidths_;
75  vector<float> bottomWidths_;
76  vector<float> lengths_;
77  vector<float> thicknesses_;
78  string infileName_;
79  string outfileName_;
81 };
82 
84  : dtGeometryToken_{esConsumes<DTGeometry, MuonGeometryRecord>(edm::ESInputTag{})},
85  infileName_(iConfig.getUntrackedParameter<string>("infileName", "cmsGeom10.root")),
86  outfileName_(iConfig.getUntrackedParameter<string>("outfileName", "validateDTGeometry.root")),
87  tolerance_(iConfig.getUntrackedParameter<int>("tolerance", 6)) {
89  outFile_ = new TFile(outfileName_.c_str(), "RECREATE");
90 }
91 
94 
95  if (dtGeometry_.isValid()) {
96  LogVerbatim("DTGeometry") << "Validating DT chamber geometry";
98 
99  LogVerbatim("DTGeometry") << "Validating DT layer geometry";
101  } else
102  LogVerbatim("DTGeometry") << "Invalid DT geometry";
103 }
104 
106  clearData();
107 
108  for (auto const& it : dtGeometry_->chambers()) {
109  DTChamberId chId = it->id();
110  GlobalPoint gp = it->surface().toGlobal(LocalPoint(0.0, 0.0, 0.0));
111 
112  const TGeoMatrix* matrix = fwGeometry_.getMatrix(chId.rawId());
113 
114  if (!matrix) {
115  LogVerbatim("DTGeometry") << "Failed to get matrix of DT chamber with detid: " << chId.rawId();
116  continue;
117  }
118 
119  compareTransform(gp, matrix);
120 
121  auto const& shape = fwGeometry_.getShapePars(chId.rawId());
122 
123  if (!shape) {
124  LogVerbatim("DTGeometry") << "Failed to get shape of DT chamber with detid: " << chId.rawId();
125  continue;
126  }
127 
128  compareShape(it, shape);
129  }
130 
131  makeHistograms("DT Chamber");
132 }
133 
135  clearData();
136 
137  vector<float> wire_positions;
138 
139  for (auto const& it : dtGeometry_->layers()) {
140  DTLayerId layerId = it->id();
141  GlobalPoint gp = it->surface().toGlobal(LocalPoint(0.0, 0.0, 0.0));
142 
143  const TGeoMatrix* matrix = fwGeometry_.getMatrix(layerId.rawId());
144 
145  if (!matrix) {
146  LogVerbatim("DTGeometry") << "Failed to get matrix of DT layer with detid: " << layerId.rawId();
147  continue;
148  }
149 
150  compareTransform(gp, matrix);
151 
152  auto const& shape = fwGeometry_.getShapePars(layerId.rawId());
153 
154  if (!shape) {
155  LogVerbatim("DTGeometry") << "Failed to get shape of DT layer with detid: " << layerId.rawId();
156  continue;
157  }
158 
159  compareShape(it, shape);
160 
161  auto const& parameters = fwGeometry_.getParameters(layerId.rawId());
162 
163  if (parameters == nullptr) {
164  LogVerbatim("DTGeometry") << "Parameters empty for DT layer with detid: " << layerId.rawId();
165  continue;
166  }
167 
168  float width = it->surface().bounds().width();
169  assert(width == parameters[6]);
170 
171  float thickness = it->surface().bounds().thickness();
172  assert(thickness == parameters[7]);
173 
174  float length = it->surface().bounds().length();
175  assert(length == parameters[8]);
176 
177  int firstChannel = it->specificTopology().firstChannel();
178  assert(firstChannel == parameters[3]);
179 
180  int lastChannel = it->specificTopology().lastChannel();
181  int nChannels = parameters[5];
182  assert(nChannels == (lastChannel - firstChannel) + 1);
183 
184  for (int wireN = firstChannel; wireN - lastChannel <= 0; ++wireN) {
185  float localX1 = it->specificTopology().wirePosition(wireN);
186  float localX2 = (wireN - (firstChannel - 1) - 0.5f) * parameters[0] - nChannels / 2.0f * parameters[0];
187  wire_positions.emplace_back(getDiff(localX1, localX2));
188  }
189  }
190 
191  makeHistogram("DT Layer Wire localX", wire_positions);
192  makeHistograms("DT Layer");
193 }
194 
196  double local[3] = {0.0, 0.0, 0.0};
197  double global[3];
198 
199  matrix->LocalToMaster(local, global);
200 
201  float distance = getDistance(GlobalPoint(global[0], global[1], global[2]), gp);
202  globalDistances_.push_back(distance);
203 }
204 
205 void DTGeometryValidate::compareShape(const GeomDet* det, const float* shape) {
206  float shapeTopWidth;
207  float shapeBottomWidth;
208  float shapeLength;
209  float shapeThickness;
210 
211  if (shape[0] == 1) {
212  shapeTopWidth = shape[2];
213  shapeBottomWidth = shape[1];
214  shapeLength = shape[4];
215  shapeThickness = shape[3];
216  } else if (shape[0] == 2) {
217  shapeTopWidth = shape[1];
218  shapeBottomWidth = shape[1];
219  shapeLength = shape[2];
220  shapeThickness = shape[3];
221  } else {
222  LogVerbatim("DTGeometry") << "Failed to get box or trapezoid from shape";
223  return;
224  }
225 
226  float topWidth, bottomWidth;
227  float length, thickness;
228 
229  const Bounds* bounds = &(det->surface().bounds());
230 
231  if (const TrapezoidalPlaneBounds* tpbs = dynamic_cast<const TrapezoidalPlaneBounds*>(bounds)) {
232  array<const float, 4> const& ps = tpbs->parameters();
233 
234  assert(ps.size() == 4);
235 
236  bottomWidth = ps[0];
237  topWidth = ps[1];
238  thickness = ps[2];
239  length = ps[3];
240  } else if ((dynamic_cast<const RectangularPlaneBounds*>(bounds))) {
241  length = det->surface().bounds().length() * 0.5;
242  topWidth = det->surface().bounds().width() * 0.5;
243  bottomWidth = topWidth;
244  thickness = det->surface().bounds().thickness() * 0.5;
245  } else {
246  LogVerbatim("DTGeometry") << "Failed to get bounds";
247  return;
248  }
249  topWidths_.push_back(fabs(shapeTopWidth - topWidth));
250  bottomWidths_.push_back(fabs(shapeBottomWidth - bottomWidth));
251  lengths_.push_back(fabs(shapeLength - length));
252  thicknesses_.push_back(fabs(shapeThickness - thickness));
253 }
254 
256  return sqrt((p1.x() - p2.x()) * (p1.x() - p2.x()) + (p1.y() - p2.y()) * (p1.y() - p2.y()) +
257  (p1.z() - p2.z()) * (p1.z() - p2.z()));
258 }
259 
260 float DTGeometryValidate::getDiff(const float val1, const float val2) {
261  if (almost_equal(val1, val2, tolerance_))
262  return 0.0f;
263  else
264  return (val1 - val2);
265 }
266 
268  outFile_->cd();
269 
270  string d(detector);
271 
272  string gdn = d + ": distance between points in global coordinates";
274 
275  string twn = d + ": absolute difference between top widths (along X)";
277 
278  string bwn = d + ": absolute difference between bottom widths (along X)";
280 
281  string ln = d + ": absolute difference between lengths (along Y)";
282  makeHistogram(ln, lengths_);
283 
284  string tn = d + ": absolute difference between thicknesses (along Z)";
286 }
287 
288 void DTGeometryValidate::makeHistogram(const string& name, vector<float>& data) {
289  if (data.empty())
290  return;
291 
292  const auto [minE, maxE] = minmax_element(begin(data), end(data));
293 
294  TH1D hist(name.c_str(), name.c_str(), 100, *minE * (1 + 0.10), *maxE * (1 + 0.10));
295 
296  for (auto const& it : data)
297  hist.Fill(it);
298 
299  hist.GetXaxis()->SetTitle("[cm]");
300  hist.Write();
301 }
302 
304 
306  LogVerbatim("DTGeometry") << "Done.";
307  LogVerbatim("DTGeometry") << "Results written to " << outfileName_;
308  outFile_->Close();
309 }
310 
type
Definition: HCALResponse.h:21
enable_if<!numeric_limits< T >::is_integer, bool >::type almost_equal(T x, T y, int ulp)
virtual float length() const =0
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
~DTGeometryValidate() override
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:474
float getDistance(const GlobalPoint &, const GlobalPoint &)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
DTGeometryValidate(const ParameterSet &)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T y() const
Definition: PV3DBase.h:60
float getDiff(const float, const float)
const Bounds & bounds() const
Definition: Surface.h:89
void makeHistograms(const char *)
void endJob() override
const TGeoMatrix * getMatrix(unsigned int id) const
Definition: FWGeometry.cc:221
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
void compareShape(const GeomDet *, const float *)
void beginJob()
Definition: Breakpoints.cc:14
virtual float width() const =0
const float * getShapePars(unsigned int id) const
Definition: FWGeometry.cc:485
void compareTransform(const GlobalPoint &, const TGeoMatrix *)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeometryToken_
vector< float > globalDistances_
edm::ESHandle< DTGeometry > dtGeometry_
void loadMap(const char *fileName)
Definition: FWGeometry.cc:99
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
vector< float > topWidths_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
#define end
Definition: vmac.h:39
T min(T a, T b)
Definition: MathUtil.h:58
double p2[4]
Definition: TauolaWrapper.h:90
d
Definition: ztail.py:151
vector< float > lengths_
vector< float > thicknesses_
void beginJob() override
virtual float thickness() const =0
const std::vector< const DTLayer * > & layers() const
Return a vector of all SuperLayer.
Definition: DTGeometry.cc:88
#define begin
Definition: vmac.h:32
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
double p1[4]
Definition: TauolaWrapper.h:89
vector< float > bottomWidths_
Definition: Bounds.h:20
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:141
bool isValid() const
Definition: ESHandle.h:44
void makeHistogram(const string &, vector< float > &)
long double T
T x() const
Definition: PV3DBase.h:59
Definition: event.py:1
void analyze(const edm::Event &, const edm::EventSetup &) override