CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ValidateGeometry.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // $Id: ValidateGeometry.cc,v 1.35 2011/10/18 09:49:30 yana Exp $
4 //
5 
13 
15 
20 
23 
27 
33 
35 
39 
43 
49 
52 
53 #include <string>
54 #include <iostream>
55 
56 #include <TEveGeoNode.h>
57 #include <TGeoVolume.h>
58 #include <TGeoBBox.h>
59 #include <TGeoArb8.h>
60 #include <TFile.h>
61 #include <TH1.h>
62 
64 {
65 public:
66  explicit ValidateGeometry(const edm::ParameterSet&);
68 
69 private:
70  virtual void beginJob();
71  virtual void analyze(const edm::Event&, const edm::EventSetup&);
72  virtual void endJob();
73 
74  void validateRPCGeometry(const int regionNumber,
75  const char* regionName);
76 
79 
80  void validateCSChamberGeometry(const int endcap,
81  const char* detname);
82 
83  void validateCSCLayerGeometry(const int endcap,
84  const char* detname);
85 
86  void validateCaloGeometry(DetId::Detector detector, int subdetector,
87  const char* detname);
88 
90  const char* detname);
91 
93  const char* detname);
94 
96  const char* detname);
97 
99  const char* detname);
100 
101  void compareTransform(const GlobalPoint& point, const TGeoMatrix* matrix);
102 
103  void compareShape(const GeomDet* det, const float* shape);
104 
105  double getDistance(const GlobalPoint& point1, const GlobalPoint& point2);
106 
107  void makeHistograms(const char* detector);
108  void makeHistogram(const std::string& name, std::vector<double>& data);
109 
110  std::string infileName_;
111  std::string outfileName_;
112 
118 
120 
121  TFile* outFile_;
122 
123  std::vector<double> globalDistances_;
124  std::vector<double> topWidths_;
125  std::vector<double> bottomWidths_;
126  std::vector<double> lengths_;
127  std::vector<double> thicknesses_;
128 
129  void clearData()
130  {
131  globalDistances_.clear();
132  topWidths_.clear();
133  bottomWidths_.clear();
134  lengths_.clear();
135  thicknesses_.clear();
136  }
137 
138  bool dataEmpty()
139  {
140  return (globalDistances_.empty() &&
141  topWidths_.empty() &&
142  bottomWidths_.empty() &&
143  lengths_.empty() &&
144  thicknesses_.empty());
145  }
146 };
147 
148 
150  : infileName_(iConfig.getUntrackedParameter<std::string>("infileName")),
151  outfileName_(iConfig.getUntrackedParameter<std::string>("outfileName"))
152 {
154 
155  outFile_ = new TFile(outfileName_.c_str(), "RECREATE");
156 }
157 
158 
160 {}
161 
162 
163 void
165 {
166  eventSetup.get<MuonGeometryRecord>().get(rpcGeometry_);
167 
168  if ( rpcGeometry_.isValid() )
169  {
170  std::cout<<"Validating RPC -z endcap geometry"<<std::endl;
171  validateRPCGeometry(-1, "RPC -z endcap");
172 
173  std::cout<<"Validating RPC +z endcap geometry"<<std::endl;
174  validateRPCGeometry(+1, "RPC +z endcap");
175 
176  std::cout<<"Validating RPC barrel geometry"<<std::endl;
177  validateRPCGeometry(0, "RPC barrel");
178  }
179  else
180  fwLog(fwlog::kWarning)<<"Invalid RPC geometry"<<std::endl;
181 
182 
183  eventSetup.get<MuonGeometryRecord>().get(dtGeometry_);
184 
185  if ( dtGeometry_.isValid() )
186  {
187  std::cout<<"Validating DT chamber geometry"<<std::endl;
189 
190  std::cout<<"Validating DT layer geometry"<<std::endl;
192  }
193  else
194  fwLog(fwlog::kWarning)<<"Invalid DT geometry"<<std::endl;
195 
196 
197  eventSetup.get<MuonGeometryRecord>().get(cscGeometry_);
198 
199  if ( cscGeometry_.isValid() )
200  {
201  std::cout<<"Validating CSC -z geometry"<<std::endl;
202  validateCSChamberGeometry(-1, "CSC chamber -z endcap");
203 
204  std::cout<<"Validating CSC +z geometry"<<std::endl;
205  validateCSChamberGeometry(+1, "CSC chamber +z endcap");
206 
207  std::cout<<"Validating CSC layer -z geometry"<<std::endl;
208  validateCSCLayerGeometry(-1, "CSC layer -z endcap");
209 
210  std::cout<<"Validating CSC layer +z geometry"<<std::endl;
211  validateCSCLayerGeometry(+1, "CSC layer +z endcap");
212  }
213  else
214  fwLog(fwlog::kWarning)<<"Invalid CSC geometry"<<std::endl;
215 
216 
218 
219  if ( trackerGeometry_.isValid() )
220  {
221  std::cout<<"Validating TIB geometry and topology"<<std::endl;
222  validateTrackerGeometry(trackerGeometry_->detsTIB(), "TIB");
223  validateStripTopology(trackerGeometry_->detsTIB(), "TIB");
224 
225  std::cout<<"Validating TOB geometry and topology"<<std::endl;
226  validateTrackerGeometry(trackerGeometry_->detsTOB(), "TOB");
227  validateStripTopology(trackerGeometry_->detsTOB(), "TOB");
228 
229  std::cout<<"Validating TEC geometry and topology"<<std::endl;
230  validateTrackerGeometry(trackerGeometry_->detsTEC(), "TEC");
231  validateStripTopology(trackerGeometry_->detsTEC(), "TEC");
232 
233  std::cout<<"Validating TID geometry and topology"<<std::endl;
234  validateTrackerGeometry(trackerGeometry_->detsTID(), "TID");
235  validateStripTopology(trackerGeometry_->detsTID(), "TID");
236 
237  std::cout<<"Validating PXB geometry and topology"<<std::endl;
238  validateTrackerGeometry(trackerGeometry_->detsPXB(), "PXB");
239  validatePixelTopology(trackerGeometry_->detsPXB(), "PXB");
240 
241  std::cout<<"Validating PXF geometry and topology"<<std::endl;
242  validateTrackerGeometry(trackerGeometry_->detsPXF(), "PXF");
243  validatePixelTopology(trackerGeometry_->detsPXF(), "PXF");
244  }
245  else
246  fwLog(fwlog::kWarning)<<"Invalid Tracker geometry"<<std::endl;
247 
248  eventSetup.get<CaloGeometryRecord>().get(caloGeometry_);
249 
250 
251  if ( caloGeometry_.isValid() )
252  {
253  std::cout<<"Validating EB geometry"<<std::endl;
255 
256  std::cout<<"Validating EE geometry"<<std::endl;
258 
259  std::cout<<"Validating ES geometry"<<std::endl;
261 
262  std::cout<<"Validating HB geometry"<<std::endl;
264 
265  std::cout<<"Validating HE geometry"<<std::endl;
267 
268  std::cout<<"Validating HO geometry"<<std::endl;
270 
271  std::cout<<"Validating HF geometry"<<std::endl;
273 
274  std::cout<<"Validating Castor geometry"<<std::endl;
276 
277  std::cout<<"Validating ZDC geometry"<<std::endl;
279  }
280  else
281  fwLog(fwlog::kWarning)<<"Invalid Calo geometry"<<std::endl;
282 
283 }
284 
285 
286 void
287 ValidateGeometry::validateRPCGeometry(const int regionNumber, const char* regionName)
288 {
289  clearData();
290 
291  std::vector<double> centers;
292 
293  std::vector<RPCRoll*> rolls = rpcGeometry_->rolls();
294 
295  for ( std::vector<RPCRoll*>::const_iterator it = rolls.begin(),
296  itEnd = rolls.end();
297  it != itEnd; ++it )
298  {
299  const RPCRoll* roll = *it;
300 
301  if ( roll )
302  {
303  RPCDetId rpcDetId = roll->id();
304 
305  if ( rpcDetId.region() == regionNumber )
306  {
307  const GeomDetUnit* det = rpcGeometry_->idToDetUnit(rpcDetId);
308  GlobalPoint gp = det->surface().toGlobal(LocalPoint(0.0, 0.0, 0.0));
309 
310  const TGeoMatrix* matrix = fwGeometry_.getMatrix(rpcDetId.rawId());
311 
312  if ( ! matrix )
313  {
314  std::cout<<"Failed to get matrix of RPC with detid: "
315  << rpcDetId.rawId() <<std::endl;
316  continue;
317  }
318 
319  compareTransform(gp, matrix);
320 
321  const float* shape = fwGeometry_.getShapePars(rpcDetId.rawId());
322 
323  if ( ! shape )
324  {
325  std::cout<<"Failed to get shape of RPC with detid: "
326  << rpcDetId.rawId() <<std::endl;
327  continue;
328  }
329 
330  compareShape(det, shape);
331 
332  const float* parameters = fwGeometry_.getParameters(rpcDetId.rawId());
333 
334  if ( parameters == 0 )
335  {
336  std::cout<<"Parameters empty for RPC with detid: "
337  << rpcDetId.rawId() <<std::endl;
338  continue;
339  }
340 
341  // Yes, I know that below I'm comparing the equivalence
342  // of floating point numbers
343 
344  int nStrips = roll->nstrips();
345  assert(nStrips == parameters[0]);
346 
347  float stripLength = roll->specificTopology().stripLength();
348  assert(stripLength == parameters[1]);
349 
350  float pitch = roll->specificTopology().pitch();
351  assert(pitch == parameters[2]);
352 
353  float offset = -0.5*nStrips*pitch;
354 
355  for ( int strip = 1; strip <= roll->nstrips(); ++strip )
356  {
357  LocalPoint centreOfStrip1 = roll->centreOfStrip(strip);
358  LocalPoint centreOfStrip2 = LocalPoint((strip-0.5)*pitch + offset, 0.0);
359 
360  centers.push_back(centreOfStrip1.x()-centreOfStrip2.x());
361  }
362  }
363  }
364  }
365 
366  std::string hn(regionName);
367  makeHistogram(hn+": centreOfStrip", centers);
368 
369 
370  makeHistograms(regionName);
371 }
372 
373 
374 void
376 {
377  clearData();
378 
379  std::vector<DTChamber*> chambers = dtGeometry_->chambers();
380 
381  for ( std::vector<DTChamber*>::const_iterator it = chambers.begin(),
382  itEnd = chambers.end();
383  it != itEnd; ++it)
384  {
385  const DTChamber* chamber = *it;
386 
387  if ( chamber )
388  {
389  DTChamberId chId = chamber->id();
390  GlobalPoint gp = chamber->surface().toGlobal(LocalPoint(0.0, 0.0, 0.0));
391 
392  const TGeoMatrix* matrix = fwGeometry_.getMatrix(chId.rawId());
393 
394  if ( ! matrix )
395  {
396  std::cout<<"Failed to get matrix of DT chamber with detid: "
397  << chId.rawId() <<std::endl;
398  continue;
399  }
400 
401  compareTransform(gp, matrix);
402 
403  const float* shape = fwGeometry_.getShapePars(chId.rawId());
404 
405  if ( ! shape )
406  {
407  std::cout<<"Failed to get shape of DT chamber with detid: "
408  << chId.rawId() <<std::endl;
409  continue;
410  }
411 
412  compareShape(chamber, shape);
413  }
414  }
415 
416  makeHistograms("DT chamber");
417 }
418 
419 
420 void
422 {
423  clearData();
424 
425  std::vector<double> wire_positions;
426 
427  std::vector<DTLayer*> layers = dtGeometry_->layers();
428 
429  for ( std::vector<DTLayer*>::const_iterator it = layers.begin(),
430  itEnd = layers.end();
431  it != itEnd; ++it)
432  {
433  const DTLayer* layer = *it;
434 
435  if ( layer )
436  {
437  DTLayerId layerId = layer->id();
438  GlobalPoint gp = layer->surface().toGlobal(LocalPoint(0.0, 0.0, 0.0));
439 
440  const TGeoMatrix* matrix = fwGeometry_.getMatrix(layerId.rawId());
441 
442  if ( ! matrix )
443  {
444  std::cout<<"Failed to get matrix of DT layer with detid: "
445  << layerId.rawId() <<std::endl;
446  continue;
447  }
448 
449  compareTransform(gp, matrix);
450 
451  const float* shape = fwGeometry_.getShapePars(layerId.rawId());
452 
453  if ( ! shape )
454  {
455  std::cout<<"Failed to get shape of DT layer with detid: "
456  << layerId.rawId() <<std::endl;
457  continue;
458  }
459 
460  compareShape(layer, shape);
461 
462 
463  const float* parameters = fwGeometry_.getParameters(layerId.rawId());
464 
465  if ( parameters == 0 )
466  {
467  std::cout<<"Parameters empty for DT layer with detid: "
468  << layerId.rawId() <<std::endl;
469  continue;
470  }
471 
472  float width = layer->surface().bounds().width();
473  assert(width == parameters[6]);
474 
475  float thickness = layer->surface().bounds().thickness();
476  assert(thickness == parameters[7]);
477 
478  float length = layer->surface().bounds().length();
479  assert(length == parameters[8]);
480 
481  int firstChannel = layer->specificTopology().firstChannel();
482  assert(firstChannel == parameters[3]);
483 
484  int lastChannel = layer->specificTopology().lastChannel();
485  int nChannels = parameters[5];
486  assert(nChannels == (lastChannel-firstChannel)+1);
487 
488  for ( int wireN = firstChannel; wireN - lastChannel <= 0; ++wireN )
489  {
490  double localX1 = layer->specificTopology().wirePosition(wireN);
491  double localX2 = (wireN -(firstChannel-1)-0.5)*parameters[0] - nChannels/2.0*parameters[0];
492 
493  wire_positions.push_back(localX1-localX2);
494 
495  //std::cout<<"wireN, localXpos: "<< wireN <<" "<< localX1 <<" "<< localX2 <<std::endl;
496 
497  }
498  }
499  }
500 
501  makeHistogram("DT layer wire localX", wire_positions);
502 
503  makeHistograms("DT layer");
504 }
505 
506 
507 void
509 {
510  clearData();
511 
512  std::vector<CSCChamber *> chambers = cscGeometry_->chambers();
513 
514  for ( std::vector<CSCChamber*>::const_iterator it = chambers.begin(),
515  itEnd = chambers.end();
516  it != itEnd; ++it )
517  {
518  const CSCChamber* chamber = *it;
519 
520  if ( chamber && chamber->id().endcap() == endcap )
521  {
522  DetId detId = chamber->geographicalId();
523  GlobalPoint gp = chamber->surface().toGlobal(LocalPoint(0.0,0.0,0.0));
524 
525  const TGeoMatrix* matrix = fwGeometry_.getMatrix(detId.rawId());
526 
527  if ( ! matrix )
528  {
529  std::cout<<"Failed to get matrix of CSC chamber with detid: "
530  << detId.rawId() <<std::endl;
531  continue;
532  }
533 
534  compareTransform(gp, matrix);
535 
536  const float* shape = fwGeometry_.getShapePars(detId.rawId());
537 
538  if ( ! shape )
539  {
540  std::cout<<"Failed to get shape of CSC chamber with detid: "
541  << detId.rawId() <<std::endl;
542  continue;
543  }
544 
545  compareShape(chamber, shape);
546  }
547  }
548 
549  makeHistograms(detname);
550 }
551 
552 void
553 ValidateGeometry::validateCSCLayerGeometry(const int endcap, const char* detname)
554 {
555  clearData();
556  std::vector<double> strip_positions;
557  std::vector<double> wire_positions;
558 
559  std::vector<double> me11_wiresLocal;
560  std::vector<double> me12_wiresLocal;
561  std::vector<double> me13_wiresLocal;
562  std::vector<double> me14_wiresLocal;
563  std::vector<double> me21_wiresLocal;
564  std::vector<double> me22_wiresLocal;
565  std::vector<double> me31_wiresLocal;
566  std::vector<double> me32_wiresLocal;
567  std::vector<double> me41_wiresLocal;
568  std::vector<double> me42_wiresLocal;
569 
570  std::vector<CSCLayer*> layers = cscGeometry_->layers();
571 
572  for ( std::vector<CSCLayer*>::const_iterator it = layers.begin(),
573  itEnd = layers.end();
574  it != itEnd; ++it )
575  {
576  const CSCLayer* layer = *it;
577 
578  if ( layer && layer->id().endcap() == endcap )
579  {
580  DetId detId = layer->geographicalId();
581  GlobalPoint gp = layer->surface().toGlobal(LocalPoint(0.0,0.0,0.0));
582 
583  const TGeoMatrix* matrix = fwGeometry_.getMatrix(detId.rawId());
584 
585  if ( ! matrix )
586  {
587  std::cout<<"Failed to get matrix of CSC layer with detid: "
588  << detId.rawId() <<std::endl;
589  continue;
590  }
591 
592  compareTransform(gp, matrix);
593 
594  const float* shape = fwGeometry_.getShapePars(detId.rawId());
595 
596  if ( ! shape )
597  {
598  std::cout<<"Failed to get shape of CSC layer with detid: "
599  << detId.rawId() <<std::endl;
600  continue;
601  }
602 
603  compareShape(layer, shape);
604 
605  double length;
606 
607  if ( shape[0] == 1 )
608  {
609  length = shape[4];
610  }
611 
612  else
613  {
614  std::cout<<"Failed to get trapezoid from shape for CSC layer with detid: "
615  << detId.rawId() <<std::endl;
616  continue;
617  }
618 
619  const float* parameters = fwGeometry_.getParameters(detId.rawId());
620 
621  if ( parameters == 0 )
622  {
623  std::cout<<"Parameters empty for CSC layer with detid: "
624  << detId.rawId() <<std::endl;
625  continue;
626  }
627 
628  int yAxisOrientation = layer->geometry()->topology()->yAxisOrientation();
629  assert(yAxisOrientation == parameters[0]);
630 
631  float centreToIntersection = layer->geometry()->topology()->centreToIntersection();
632  assert(centreToIntersection == parameters[1]);
633 
634  float yCentre = layer->geometry()->topology()->yCentreOfStripPlane();
635  assert(yCentre == parameters[2]);
636 
637  float phiOfOneEdge = layer->geometry()->topology()->phiOfOneEdge();
638  assert(phiOfOneEdge == parameters[3]);
639 
640  float stripOffset = layer->geometry()->topology()->stripOffset();
641  assert(stripOffset == parameters[4]);
642 
643  float angularWidth = layer->geometry()->topology()->angularWidth();
644  assert(angularWidth == parameters[5]);
645 
646  for ( int nStrip = 1; nStrip <= layer->geometry()->numberOfStrips();
647  ++nStrip )
648  {
649  float xOfStrip1 = layer->geometry()->xOfStrip(nStrip);
650 
651  double stripAngle = phiOfOneEdge + yAxisOrientation*(nStrip-(0.5-stripOffset))*angularWidth;
652  double xOfStrip2 = yAxisOrientation*(centreToIntersection-yCentre)*tan(stripAngle);
653 
654  strip_positions.push_back(xOfStrip1-xOfStrip2);
655  }
656 
657  int station = layer->id().station();
658  int ring = layer->id().ring();
659 
660  double wireSpacingInGroup = layer->geometry()->wireTopology()->wireSpacing();
661  assert(wireSpacingInGroup == parameters[6]);
662 
663  double wireSpacing = 0.0;
664  // we calculate an average wire group
665  // spacing from radialExtentOfTheWirePlane / numOfWireGroups
666 
667  double extentOfWirePlane = 0.0;
668 
669  if ( ring == 2 )
670  {
671  if ( station == 1 )
672  extentOfWirePlane = 174.81; //wireSpacing = 174.81/64;
673  else
674  extentOfWirePlane = 323.38; //wireSpacing = 323.38/64;
675  }
676  else if ( station == 1 && (ring == 1 || ring == 4))
677  extentOfWirePlane = 150.5; //wireSpacing = 150.5/48;
678  else if ( station == 1 && ring == 3 )
679  extentOfWirePlane = 164.47; //wireSpacing = 164.47/32;
680  else if ( station == 2 && ring == 1 )
681  extentOfWirePlane = 189.97; //wireSpacing = 189.97/112;
682  else if ( station == 3 && ring == 1 )
683  extentOfWirePlane = 170.01; //wireSpacing = 170.01/96;
684  else if ( station == 4 && ring == 1 )
685  extentOfWirePlane = 149.73; //wireSpacing = 149.73/96;
686 
687  float wireAngle = layer->geometry()->wireTopology()->wireAngle();
688  assert(wireAngle == parameters[7]);
689 
690  //float cosWireAngle = cos(wireAngle);
691 
692  /* NOTE
693  Some parameters don't seem available in a public interface
694  so have to perhaps hard-code. This may not be too bad as there
695  seems to be a lot of degeneracy.
696  */
697 
698  double alignmentPinToFirstWire;
699  double yAlignmentFrame = 3.49;
700 
701  if ( station == 1 )
702  {
703  if ( ring == 1 || ring == 4 )
704  {
705  alignmentPinToFirstWire = 1.065;
706  yAlignmentFrame = 0.0;
707  }
708 
709  else // ME12, ME 13
710  alignmentPinToFirstWire = 2.85;
711  }
712 
713  else if ( station == 4 && ring == 1 )
714  alignmentPinToFirstWire = 3.04;
715 
716  else if ( station == 3 && ring == 1 )
717  alignmentPinToFirstWire = 2.84;
718 
719  else // ME21, ME22, ME32, ME42
720  alignmentPinToFirstWire = 2.87;
721 
722  double yOfFirstWire = (yAlignmentFrame-length) + alignmentPinToFirstWire;
723 
724  int nWireGroups = layer->geometry()->numberOfWireGroups();
725  double E = extentOfWirePlane/nWireGroups;
726 
727  for ( int nWireGroup = 1; nWireGroup <= nWireGroups; ++nWireGroup )
728  {
729  LocalPoint centerOfWireGroup = layer->geometry()->localCenterOfWireGroup(nWireGroup);
730  double yOfWire1 = centerOfWireGroup.y();
731 
732  //double yOfWire2 = (-0.5 - (nWireGroups*0.5 - 1) + (nWireGroup-1))*E;
733  //yOfWire2 += 0.5*E;
734  double yOfWire2 = yOfFirstWire + ((nWireGroup-1)*E);
735  yOfWire2 += wireSpacing*0.5;
736 
737  double ydiff_local = yOfWire1 - yOfWire2;
738  wire_positions.push_back(ydiff_local);
739 
740  //GlobalPoint globalPoint = layer->surface().toGlobal(LocalPoint(0.0,yOfWire1,0.0));
741 
742  /*
743  float fwLocalPoint[3] =
744  {
745  0.0, yOfWire2, 0.0
746  };
747 
748  float fwGlobalPoint[3];
749  fwGeometry_.localToGlobal(detId.rawId(), fwLocalPoint, fwGlobalPoint);
750  double ydiff_global = globalPoint.y() - fwGlobalPoint[1];
751  */
752 
753  if ( station == 1 )
754  {
755  if ( ring == 1 )
756  {
757  me11_wiresLocal.push_back(ydiff_local);
758  }
759  else if ( ring == 2 )
760  {
761  me12_wiresLocal.push_back(ydiff_local);
762  }
763  else if ( ring == 3 )
764  {
765  me13_wiresLocal.push_back(ydiff_local);
766  }
767  else if ( ring == 4 )
768  {
769  me14_wiresLocal.push_back(ydiff_local);
770  }
771  }
772  else if ( station == 2 )
773  {
774  if ( ring == 1 )
775  {
776  me21_wiresLocal.push_back(ydiff_local);
777  }
778  else if ( ring == 2 )
779  {
780  me22_wiresLocal.push_back(ydiff_local);
781  }
782  }
783  else if ( station == 3 )
784  {
785  if ( ring == 1 )
786  {
787  me31_wiresLocal.push_back(ydiff_local);
788  }
789  else if ( ring == 2 )
790  {
791  me32_wiresLocal.push_back(ydiff_local);
792  }
793  }
794  else if ( station == 4 )
795  {
796  if ( ring == 1 )
797  {
798  me41_wiresLocal.push_back(ydiff_local);
799  }
800  else if ( ring == 2 )
801  {
802  me42_wiresLocal.push_back(ydiff_local);
803  }
804  }
805  }
806  }
807  }
808 
809  std::string hn(detname);
810  makeHistogram(hn+": xOfStrip", strip_positions);
811 
812  makeHistogram(hn+": local yOfWire", wire_positions);
813 
814  makeHistogram("ME11: local yOfWire", me11_wiresLocal);
815  makeHistogram("ME12: local yOfWire", me12_wiresLocal);
816  makeHistogram("ME13: local yOfWire", me13_wiresLocal);
817  makeHistogram("ME14: local yOfWire", me14_wiresLocal);
818  makeHistogram("ME21: local yOfWire", me21_wiresLocal);
819  makeHistogram("ME22: local yOfWire", me22_wiresLocal);
820  makeHistogram("ME31: local yOfWire", me31_wiresLocal);
821  makeHistogram("ME32: local yOfWire", me32_wiresLocal);
822  makeHistogram("ME41: local yOfWire", me41_wiresLocal);
823  makeHistogram("ME42: local yOfWire", me42_wiresLocal);
824 
825  makeHistograms(detname);
826 }
827 
828 void
830  int subdetector,
831  const char* detname)
832 {
833  clearData();
834 
836  caloGeometry_->getSubdetectorGeometry(detector, subdetector);
837 
838  const std::vector<DetId>& ids = geometry->getValidDetIds(detector, subdetector);
839 
840  for (std::vector<DetId>::const_iterator it = ids.begin(),
841  iEnd = ids.end();
842  it != iEnd; ++it)
843  {
844  unsigned int rawId = (*it).rawId();
845 
846  const float* points = fwGeometry_.getCorners(rawId);
847 
848  if ( points == 0 )
849  {
850  std::cout <<"Failed to get points of "<< detname
851  <<" element with detid: "<< rawId <<std::endl;
852  continue;
853  }
854 
855  const CaloCellGeometry* cellGeometry = geometry->getGeometry(*it);
856  const CaloCellGeometry::CornersVec& corners = cellGeometry->getCorners();
857 
858  assert(corners.size() == 8);
859 
860  for ( unsigned int i = 0, offset = 0; i < 8; ++i )
861  {
862  offset = 2*i;
863 
864  double distance = getDistance(GlobalPoint(points[i+offset], points[i+1+offset], points[i+2+offset]),
865  GlobalPoint(corners[i].x(), corners[i].y(), corners[i].z()));
866 
867  globalDistances_.push_back(distance);
868  }
869  }
870 
871  makeHistograms(detname);
872 }
873 
874 
875 void
877  const char* detname)
878 {
879  clearData();
880 
881  for ( TrackerGeometry::DetContainer::const_iterator it = dets.begin(),
882  itEnd = dets.end();
883  it != itEnd; ++it )
884  {
885  GlobalPoint gp = (trackerGeometry_->idToDet((*it)->geographicalId()))->surface().toGlobal(LocalPoint(0.0,0.0,0.0));
886  unsigned int rawId = (*it)->geographicalId().rawId();
887 
888  const TGeoMatrix* matrix = fwGeometry_.getMatrix(rawId);
889 
890  if ( ! matrix )
891  {
892  std::cout <<"Failed to get matrix of "<< detname
893  <<" element with detid: "<< rawId <<std::endl;
894  continue;
895  }
896 
897  compareTransform(gp, matrix);
898 
899  const float* shape = fwGeometry_.getShapePars(rawId);
900 
901  if ( ! shape )
902  {
903  std::cout<<"Failed to get shape of "<< detname
904  <<" element with detid: "<< rawId <<std::endl;
905  continue;
906  }
907 
908  compareShape(*it, shape);
909  }
910 
911  makeHistograms(detname);
912 }
913 
914 
915 void
917  const char* detname)
918 {
919  clearData();
920 
921  for ( TrackerGeometry::DetUnitContainer::const_iterator it = dets.begin(),
922  itEnd = dets.end();
923  it != itEnd; ++it )
924  {
925  GlobalPoint gp = (trackerGeometry_->idToDet((*it)->geographicalId()))->surface().toGlobal(LocalPoint(0.0,0.0,0.0));
926  unsigned int rawId = (*it)->geographicalId().rawId();
927 
928  const TGeoMatrix* matrix = fwGeometry_.getMatrix(rawId);
929 
930  if ( ! matrix )
931  {
932  std::cout<< "Failed to get matrix of "<< detname
933  <<" element with detid: "<< rawId <<std::endl;
934  continue;
935  }
936 
937  compareTransform(gp, matrix);
938 
939 
940  const float* shape = fwGeometry_.getShapePars(rawId);
941 
942  if ( ! shape )
943  {
944  std::cout<<"Failed to get shape of "<< detname
945  <<" element with detid: "<< rawId <<std::endl;
946  continue;
947  }
948 
949  compareShape(*it, shape);
950  }
951 
952  makeHistograms(detname);
953 }
954 
955 void
957  const char* detname)
958 {
959  std::vector<double> pixelLocalXs;
960  std::vector<double> pixelLocalYs;
961 
962  for ( TrackerGeometry::DetContainer::const_iterator it = dets.begin(),
963  itEnd = dets.end();
964  it != itEnd; ++it )
965  {
966  unsigned int rawId = (*it)->geographicalId().rawId();
967 
968  const float* parameters = fwGeometry_.getParameters(rawId);
969 
970  if ( parameters == 0 )
971  {
972  std::cout<<"Parameters empty for "<< detname <<" element with detid: "
973  << rawId <<std::endl;
974  continue;
975  }
976 
977  if ( const PixelGeomDetUnit* det =
978  dynamic_cast<const PixelGeomDetUnit*>(trackerGeometry_->idToDetUnit((*it)->geographicalId())) )
979  {
980  if ( const PixelTopology* rpt = &det->specificTopology() )
981  {
982  int nrows = rpt->nrows();
983  int ncolumns = rpt->ncolumns();
984 
985  assert(parameters[0] == nrows);
986  assert(parameters[1] == ncolumns);
987 
988  for ( int row = 1; row <= nrows; ++row )
989  {
990  for ( int column = 1; column <= ncolumns; ++column )
991  {
992  LocalPoint localPoint = rpt->localPosition(MeasurementPoint(row, column));
993 
994  pixelLocalXs.push_back(localPoint.x() - fireworks::pixelLocalX(row, nrows));
995  pixelLocalYs.push_back(localPoint.y() - fireworks::pixelLocalY(column, ncolumns));
996  }
997  }
998  }
999 
1000  else
1001  std::cout<<"No topology for "<< detname <<" "<< rawId <<std::endl;
1002  }
1003 
1004  else
1005  std::cout<<"No geomDetUnit for "<< detname <<" "<< rawId <<std::endl;
1006  }
1007 
1008  std::string hn(detname);
1009  makeHistogram(hn+" pixelLocalX", pixelLocalXs);
1010  makeHistogram(hn+" pixelLocalY", pixelLocalYs);
1011 }
1012 
1013 void
1015  const char* detname)
1016 {
1017  std::vector<double> radialStripLocalXs;
1018  std::vector<double> rectangularStripLocalXs;
1019 
1020  for ( TrackerGeometry::DetContainer::const_iterator it = dets.begin(),
1021  itEnd = dets.end();
1022  it != itEnd; ++it )
1023  {
1024  unsigned int rawId = (*it)->geographicalId().rawId();
1025 
1026  const float* parameters = fwGeometry_.getParameters(rawId);
1027 
1028  if ( parameters == 0 )
1029  {
1030  std::cout<<"Parameters empty for "<< detname <<" element with detid: "
1031  << rawId <<std::endl;
1032  continue;
1033  }
1034 
1035  if ( const StripGeomDetUnit* det =
1036  dynamic_cast<const StripGeomDetUnit*>(trackerGeometry_->idToDet((*it)->geographicalId())) )
1037  {
1038  // NOTE: why the difference in dets vs. units between these and pixels? The dynamic cast above
1039  // fails for many of the detids...
1040 
1041  const StripTopology* st = dynamic_cast<const StripTopology*>(&det->specificTopology());
1042 
1043  if ( st )
1044  {
1045  //assert(parameters[0] == 0);
1046  int nstrips = st->nstrips();
1047  assert(parameters[1] == nstrips);
1048  assert(parameters[2] == st->stripLength());
1049 
1050  if( const RadialStripTopology* rst = dynamic_cast<const RadialStripTopology*>(&(det->specificType().specificTopology())) )
1051  {
1052  assert(parameters[0] == 1);
1053  assert(parameters[3] == rst->yAxisOrientation());
1054  assert(parameters[4] == rst->originToIntersection());
1055  assert(parameters[5] == rst->phiOfOneEdge());
1056  assert(parameters[6] == rst->angularWidth());
1057 
1058  for ( uint16_t strip = 1; strip <= nstrips; ++strip )
1059  {
1060  float stripAngle1 = rst->stripAngle(strip);
1061  float stripAngle2 = parameters[5] + parameters[3]*strip*parameters[6];
1062 
1063  assert((stripAngle1-stripAngle2) == 0);
1064 
1065  LocalPoint stripPosition = st->localPosition(strip);
1066 
1067  float stripX = parameters[3]*parameters[4]*tan(stripAngle2);
1068  radialStripLocalXs.push_back(stripPosition.x()-stripX);
1069  }
1070  }
1071 
1072  else if( dynamic_cast<const RectangularStripTopology*>(&(det->specificType().specificTopology())) )
1073  {
1074  assert(parameters[0] == 2);
1075  assert(parameters[3] == st->pitch());
1076 
1077  for ( uint16_t strip = 1; strip <= nstrips; ++strip )
1078  {
1079  LocalPoint stripPosition = st->localPosition(strip);
1080  float stripX = -parameters[1]*0.5*parameters[3];
1081  stripX += strip*parameters[3];
1082  rectangularStripLocalXs.push_back(stripPosition.x()-stripX);
1083  }
1084  }
1085 
1086  else if( dynamic_cast<const TrapezoidalStripTopology*>(&(det->specificType().specificTopology())) )
1087  {
1088  assert(parameters[0] == 3);
1089  assert(parameters[3] == st->pitch());
1090  }
1091 
1092  else
1093  std::cout<<"Failed to get pitch for "<< detname <<" "<< rawId <<std::endl;
1094  }
1095 
1096  else
1097  std::cout<<"Failed cast to StripTopology for "<< detname <<" "<< rawId <<std::endl;
1098  }
1099 
1100  //else
1101  // std::cout<<"Failed cast to StripGeomDetUnit for "<< detname <<" "<< rawId <<std::endl;
1102  }
1103 
1104  std::string hn(detname);
1105  makeHistogram(hn+" radial strip localX", radialStripLocalXs);
1106  makeHistogram(hn+" rectangular strip localX", rectangularStripLocalXs);
1107 }
1108 
1109 
1110 void
1112  const TGeoMatrix* matrix)
1113 {
1114  double local[3] =
1115  {
1116  0.0, 0.0, 0.0
1117  };
1118 
1119  double global[3];
1120 
1121  matrix->LocalToMaster(local, global);
1122 
1123  double distance = getDistance(GlobalPoint(global[0], global[1], global[2]), gp);
1124  globalDistances_.push_back(distance);
1125 }
1126 
1127 void
1128 ValidateGeometry::compareShape(const GeomDet* det, const float* shape)
1129 {
1130  double shape_topWidth;
1131  double shape_bottomWidth;
1132  double shape_length;
1133  double shape_thickness;
1134 
1135  if ( shape[0] == 1 )
1136  {
1137  shape_topWidth = shape[2];
1138  shape_bottomWidth = shape[1];
1139  shape_length = shape[4];
1140  shape_thickness = shape[3];
1141  }
1142 
1143  else if ( shape[0] == 2 )
1144  {
1145  shape_topWidth = shape[1];
1146  shape_bottomWidth = shape[1];
1147  shape_length = shape[2];
1148  shape_thickness = shape[3];
1149  }
1150 
1151  else
1152  {
1153  std::cout<<"Failed to get box or trapezoid from shape"<<std::endl;
1154  return;
1155  }
1156 
1157  double topWidth, bottomWidth;
1158  double length, thickness;
1159 
1160  const Bounds* bounds = &(det->surface().bounds());
1161 
1162  if ( const TrapezoidalPlaneBounds* tpbs = dynamic_cast<const TrapezoidalPlaneBounds*>(bounds) )
1163  {
1164  std::vector<float> ps = tpbs->parameters();
1165 
1166  assert(ps.size() == 4);
1167 
1168  bottomWidth = ps[0];
1169  topWidth = ps[1];
1170  thickness = ps[2];
1171  length = ps[3];
1172  }
1173 
1174  else if ( (dynamic_cast<const RectangularPlaneBounds*>(bounds)) )
1175  {
1176  length = det->surface().bounds().length()*0.5;
1177  topWidth = det->surface().bounds().width()*0.5;
1178  bottomWidth = topWidth;
1179  thickness = det->surface().bounds().thickness()*0.5;
1180  }
1181 
1182  else
1183  {
1184  std::cout<<"Failed to get bounds"<<std::endl;
1185  return;
1186  }
1187 
1188  //assert((tgeotrap && trapezoid) || (tgeobbox && rectangle));
1189 
1190  /*
1191  std::cout<<"topWidth: "<< shape_topWidth <<" "<< topWidth <<std::endl;
1192  std::cout<<"bottomWidth: "<< shape_bottomWidth <<" "<< bottomWidth <<std::endl;
1193  std::cout<<"length: "<< shape_length <<" "<< length <<std::endl;
1194  std::cout<<"thickness: "<< shape_thickness <<" "<< thickness <<std::endl;
1195  */
1196 
1197  topWidths_.push_back(fabs(shape_topWidth - topWidth));
1198  bottomWidths_.push_back(fabs(shape_bottomWidth - bottomWidth));
1199  lengths_.push_back(fabs(shape_length - length));
1200  thicknesses_.push_back(fabs(shape_thickness - thickness));
1201 
1202  return;
1203 }
1204 
1205 
1206 
1207 
1208 double
1210 {
1211  /*
1212  std::cout<<"X: "<< p1.x() <<" "<< p2.x() <<std::endl;
1213  std::cout<<"Y: "<< p1.y() <<" "<< p2.y() <<std::endl;
1214  std::cout<<"Z: "<< p1.z() <<" "<< p2.z() <<std::endl;
1215  */
1216 
1217  return sqrt((p1.x()-p2.x())*(p1.x()-p2.x())+
1218  (p1.y()-p2.y())*(p1.y()-p2.y())+
1219  (p1.z()-p2.z())*(p1.z()-p2.z()));
1220 }
1221 
1222 
1223 void
1225 {
1226  outFile_->cd();
1227 
1228  std::string d(detector);
1229 
1230  std::string gdn = d+": distance between points in global coordinates";
1232 
1233  std::string twn = d + ": absolute difference between top widths (along X)";
1234  makeHistogram(twn, topWidths_);
1235 
1236  std::string bwn = d + ": absolute difference between bottom widths (along X)";
1238 
1239  std::string ln = d + ": absolute difference between lengths (along Y)";
1240  makeHistogram(ln, lengths_);
1241 
1242  std::string tn = d + ": absolute difference between thicknesses (along Z)";
1244 
1245  return;
1246 }
1247 
1248 
1249 void
1250 ValidateGeometry::makeHistogram(const std::string& name, std::vector<double>& data)
1251 {
1252  if ( data.empty() )
1253  return;
1254 
1255  std::vector<double>::iterator it;
1256 
1257  it = std::min_element(data.begin(), data.end());
1258  double minE = *it;
1259 
1260  it = std::max_element(data.begin(), data.end());
1261  double maxE = *it;
1262 
1263  std::vector<double>::iterator itEnd = data.end();
1264 
1265  TH1D hist(name.c_str(), name.c_str(), 100, minE*(1+0.10), maxE*(1+0.10));
1266 
1267  for ( it = data.begin(); it != itEnd; ++it )
1268  hist.Fill(*it);
1269 
1270  hist.GetXaxis()->SetTitle("[cm]");
1271  hist.Write();
1272 }
1273 
1274 void
1276 {
1277  outFile_->cd();
1278 }
1279 
1280 
1281 void
1283 {
1284  std::cout<<"Done. "<<std::endl;
1285  std::cout<<"Results written to "<< outfileName_ <<std::endl;
1286  outFile_->Close();
1287 }
1288 
1290 
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:78
virtual int nstrips() const =0
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:88
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
edm::ESHandle< CaloGeometry > caloGeometry_
LocalPoint centreOfStrip(int strip) const
Definition: RPCRoll.cc:52
edm::ESHandle< TrackerGeometry > trackerGeometry_
const CSCWireTopology * wireTopology() const
virtual float length() const =0
void strip(std::string &input, const std::string &blanks=" \n\t")
Definition: stringTools.cc:16
std::vector< double > topWidths_
void validateTrackerGeometry(const TrackerGeometry::DetContainer &dets, const char *detname)
float pixelLocalY(const double mpy, const int m_ncols)
Definition: TrackUtils.cc:198
CSCDetId id() const
Get the (concrete) DetId.
Definition: CSCChamber.h:37
const float * getParameters(unsigned int id) const
Definition: FWGeometry.cc:293
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
DTLayerId id() const
Return the DetId of this SL.
Definition: DTLayer.cc:48
CSCDetId id() const
Definition: CSCLayer.h:42
edm::ESHandle< RPCGeometry > rpcGeometry_
virtual void beginJob()
int nstrips() const
Definition: RPCRoll.cc:46
std::vector< double > bottomWidths_
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
void makeHistograms(const char *detector)
T y() const
Definition: PV3DBase.h:62
float centreToIntersection() const
static char chambers[TOTALCHAMBERS][20]
Definition: ReadPGInfo.cc:243
int numberOfStrips() const
std::vector< GeomDet * > DetContainer
double wireSpacing() const
const TGeoMatrix * getMatrix(unsigned int id) const
Definition: FWGeometry.cc:166
const StripTopology & specificTopology() const
Definition: RPCRoll.cc:107
int numberOfWireGroups() const
double double double z
void makeHistogram(const std::string &name, std::vector< double > &data)
virtual const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const
Get a list of valid detector ids (for the given subdetector)
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
int firstChannel() const
Returns the wire number of the first wire.
Definition: DTTopology.h:80
float xOfStrip(int strip, float y=0.) const
void validateDTChamberGeometry()
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
int endcap() const
Definition: CSCDetId.h:95
virtual float thickness() const =0
virtual void endJob()
const float * getShapePars(unsigned int id) const
Definition: FWGeometry.cc:309
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:82
const DTTopology & specificTopology() const
Definition: DTLayer.cc:44
virtual float stripLength() const =0
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
RPCDetId id() const
Definition: RPCRoll.cc:24
float wireAngle() const
void compareShape(const GeomDet *det, const float *shape)
void loadMap(const char *fileName)
Definition: FWGeometry.cc:64
std::vector< GeomDetUnit * > DetUnitContainer
T sqrt(T t)
Definition: SSEVec.h:46
void validateCaloGeometry(DetId::Detector detector, int subdetector, const char *detname)
DTChamberId id() const
Return the DTChamberId of this chamber.
Definition: DTChamber.cc:35
T z() const
Definition: PV3DBase.h:63
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
edm::ESHandle< DTGeometry > dtGeometry_
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:72
void validateRPCGeometry(const int regionNumber, const char *regionName)
void validatePixelTopology(const TrackerGeometry::DetContainer &dets, const char *detname)
unsigned int offset(bool)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
void validateCSCLayerGeometry(const int endcap, const char *detname)
static const int SubdetectorId
double p2[4]
Definition: TauolaWrapper.h:90
edm::ESHandle< CSCGeometry > cscGeometry_
LocalPoint localCenterOfWireGroup(int wireGroup) const
int ring() const
Definition: CSCDetId.h:77
Definition: DetId.h:20
const Bounds & bounds() const
Definition: BoundSurface.h:89
const CSCStripTopology * topology() const
float yCentreOfStripPlane() const
size_type size() const
Definition: EZArrayFL.h:81
virtual float stripOffset(void) const
static const int SubdetectorId
Definition: HcalZDCDetId.h:22
const float * getCorners(unsigned int id) const
Definition: FWGeometry.cc:277
double getDistance(const GlobalPoint &point1, const GlobalPoint &point2)
#define fwLog(_level_)
Definition: fwLog.h:51
Detector
Definition: DetId.h:26
const T & get() const
Definition: EventSetup.h:55
ValidateGeometry(const edm::ParameterSet &)
ESHandle< TrackerGeometry > geometry
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double p1[4]
Definition: TauolaWrapper.h:89
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
std::vector< double > thicknesses_
virtual float pitch() const =0
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
virtual void analyze(const edm::Event &, const edm::EventSetup &)
virtual LocalPoint localPosition(float strip) const =0
int station() const
Definition: CSCDetId.h:88
tuple cout
Definition: gather_cfg.py:121
float angularWidth() const
std::string infileName_
Definition: Bounds.h:18
float phiOfOneEdge() const
Definition: DDAxes.h:10
void validateCSChamberGeometry(const int endcap, const char *detname)
void validateStripTopology(const TrackerGeometry::DetContainer &dets, const char *detname)
bool isValid() const
Definition: ESHandle.h:37
std::vector< double > lengths_
T x() const
Definition: PV3DBase.h:61
virtual float width() const =0
const CSCLayerGeometry * geometry() const
Definition: CSCLayer.h:47
float pixelLocalX(const double mpx, const int m_nrows)
Definition: TrackUtils.cc:165
std::string outfileName_
void compareTransform(const GlobalPoint &point, const TGeoMatrix *matrix)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
virtual const CornersVec & getCorners() const =0
Returns the corner points of this cell&#39;s volume.
std::vector< double > globalDistances_
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:67