CMS 3D CMS Logo

TrackerAlignment_PayloadInspector.cc
Go to the documentation of this file.
1 
10 
14 
15 // the data format of the condition to be inspected
20 
21 //#include "CondFormats/Alignment/interface/Definitions.h"
22 #include "CLHEP/Vector/RotationInterfaces.h"
24 
25 // needed for the tracker map
27 
28 // needed for mapping
31 
32 #include <memory>
33 #include <sstream>
34 #include <iostream>
35 
36 // include ROOT
37 #include "TH2F.h"
38 #include "TLegend.h"
39 #include "TCanvas.h"
40 #include "TLine.h"
41 #include "TStyle.h"
42 #include "TLatex.h"
43 #include "TPave.h"
44 #include "TPaveStats.h"
45 
46 namespace {
47 
48  // M.M. 2017/09/29
49  // Hardcoded Tracker Global Position Record
50  // Without accessing the ES, it is not possible to access to the GPR with the PI technology,
51  // so this needs to be hardcoded.
52  // Anyway it is not likely to change until a new Tracker is installed.
53  // Details at:
54  // - https://indico.cern.ch/event/238026/contributions/513928/attachments/400000/556192/mm_TkAlMeeting_28_03_2013.pdf
55  // - https://twiki.cern.ch/twiki/bin/view/CMS/TkAlignmentPixelPosition
56 
57  static const std::map<AlignmentPI::coordinate,float> hardcodeGPR =
58  {{AlignmentPI::t_x,-9.00e-02},
59  {AlignmentPI::t_y,-1.10e-01},
60  {AlignmentPI::t_z,-1.70e-01}};
61 
62  //*******************************************//
63  // Size of the movement over all partitions
64  //******************************************//
65 
66  template<AlignmentPI::coordinate coord> class TrackerAlignmentCompare : public cond::payloadInspector::PlotImage<Alignments> {
67  public:
68  TrackerAlignmentCompare() : cond::payloadInspector::PlotImage<Alignments>( "comparison of "+AlignmentPI::getStringFromCoordinate(coord)+" coordinate between two geometries" ){
69  setSingleIov( false );
70  }
71 
72  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
73 
74  std::vector<std::tuple<cond::Time_t,cond::Hash> > sorted_iovs = iovs;
75 
76  // make absolute sure the IOVs are sortd by since
77  std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const &t1, auto const &t2) {
78  return std::get<0>(t1) < std::get<0>(t2);
79  });
80 
81  auto firstiov = sorted_iovs.front();
82  auto lastiov = sorted_iovs.back();
83 
84  std::shared_ptr<Alignments> last_payload = fetchPayload( std::get<1>(lastiov) );
85  std::shared_ptr<Alignments> first_payload = fetchPayload( std::get<1>(firstiov) );
86 
87  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
88  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
89 
90  std::vector<AlignTransform> ref_ali = first_payload->m_align;
91  std::vector<AlignTransform> target_ali = last_payload->m_align;
92 
93  TCanvas canvas("Alignment Comparison","Alignment Comparison",1200,1200);
94 
95  if(ref_ali.size() != target_ali.size()){
96  edm::LogError("TrackerAlignment_PayloadInspector") << "the size of the reference alignment (" << ref_ali.size() << ") is different from the one of the target ("<< target_ali.size() << ")! You are probably trying to compare different underlying geometries. Exiting";
97  return false;
98  }
99 
100  // check that the geomtery is a tracker one
101  const char * path_toTopologyXML = (ref_ali.size()==AlignmentPI::phase0size) ? "Geometry/TrackerCommonData/data/trackerParameters.xml" : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
103 
104  for (const auto &ali : ref_ali){
105  auto mydetid = ali.rawId();
106  if(DetId(mydetid).det() != DetId::Tracker){
107  edm::LogWarning("TrackerAlignment_PayloadInspector") << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() <<" ("<< DetId(mydetid).det()<<") is different from "<<DetId::Tracker<<" (is DoubleSide: "<< tTopo.tidIsDoubleSide(mydetid)<<"); subdetId "<< DetId(mydetid).subdetId() <<" - terminating ";
108  return false;
109  }
110  }
111 
112  int counter=0;
113  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
114  std::string unit = (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z ) ? "[#mum]" : "[mrad]";
115 
116  //std::unique_ptr<TH1F> compare = std::unique_ptr<TH1F>(new TH1F("comparison",Form("Comparison of %s;DetId index; #Delta%s %s",s_coord.c_str(),s_coord.c_str(),unit.c_str()),ref_ali.size(),-0.5,ref_ali.size()-0.5));
117  std::unique_ptr<TH1F> compare = std::unique_ptr<TH1F>(new TH1F("comparison",Form(";Detector Id index; #Delta%s %s",s_coord.c_str(),unit.c_str()),ref_ali.size(),-0.5,ref_ali.size()-0.5));
118 
119  std::vector<int> boundaries;
121  for(unsigned int i=0;i<=ref_ali.size();i++){
122 
123  if(ref_ali[i].rawId() == target_ali[i].rawId()){
124 
125  counter++;
126  int subid = DetId(ref_ali[i].rawId()).subdetId();
127 
128  auto thePart = static_cast<AlignmentPI::partitions>(subid);
129  if( thePart != currentPart ){
130  currentPart=thePart;
131  boundaries.push_back(counter);
132  }
133 
134  CLHEP::HepRotation target_rot( target_ali[i].rotation() );
135  CLHEP::HepRotation ref_rot( ref_ali[i].rotation() );
136 
137  align::RotationType target_rotation( target_rot.xx(), target_rot.xy(), target_rot.xz(),
138  target_rot.yx(), target_rot.yy(), target_rot.yz(),
139  target_rot.zx(), target_rot.zy(), target_rot.zz() );
140 
141  align::RotationType ref_rotation( ref_rot.xx(), ref_rot.xy(), ref_rot.xz(),
142  ref_rot.yx(), ref_rot.yy(), ref_rot.yz(),
143  ref_rot.zx(), ref_rot.zy(), ref_rot.zz() );
144 
145  align::EulerAngles target_eulerAngles = align::toAngles(target_rotation);
146  align::EulerAngles ref_eulerAngles = align::toAngles(ref_rotation);
147 
148  switch(coord){
149  case AlignmentPI::t_x :
150  compare->SetBinContent(i+1,(target_ali[i].translation().x()-ref_ali[i].translation().x())*AlignmentPI::cmToUm);
151  break;
152  case AlignmentPI::t_y:
153  compare->SetBinContent(i+1,(target_ali[i].translation().y()-ref_ali[i].translation().y())*AlignmentPI::cmToUm);
154  break;
155  case AlignmentPI::t_z:
156  compare->SetBinContent(i+1,(target_ali[i].translation().z()-ref_ali[i].translation().z())*AlignmentPI::cmToUm);
157  break;
159  compare->SetBinContent(i+1,(target_eulerAngles[0]-ref_eulerAngles[0])*1000.);
160  break;
162  compare->SetBinContent(i+1,(target_eulerAngles[1]-ref_eulerAngles[1])*1000.);
163  break;
165  compare->SetBinContent(i+1,(target_eulerAngles[2]-ref_eulerAngles[2])*1000.);
166  break;
167  default:
168  edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate "<< coord << std::endl;
169  break;
170  } // switch on the coordinate
171  } // check on the same detID
172  } // loop on the components
173 
174  canvas.cd();
175 
176  canvas.SetLeftMargin(0.17);
177  canvas.SetRightMargin(0.05);
178  canvas.SetBottomMargin(0.15);
179  AlignmentPI::makeNicePlotStyle(compare.get(),kBlack);
180  auto max = compare->GetMaximum();
181  auto min = compare->GetMinimum();
182  auto range = std::abs(max) > std::abs(min) ? std::abs(max) : std::abs(min);
183  //auto newMax = (max > 0.) ? max*1.2 : max*0.8;
184  compare->GetYaxis()->SetRangeUser(-range*1.3,range*1.2);
185  compare->SetMarkerStyle(20);
186  compare->SetMarkerSize(0.5);
187  compare->Draw("P");
188 
189  canvas.Update();
190  canvas.cd();
191 
192  TLine l[boundaries.size()];
193  unsigned int i=0;
194  for (const auto & line : boundaries){
195  l[i] = TLine(compare->GetBinLowEdge(line),canvas.cd()->GetUymin(),compare->GetBinLowEdge(line),canvas.cd()->GetUymax());
196  l[i].SetLineWidth(1);
197  l[i].SetLineStyle(9);
198  l[i].SetLineColor(2);
199  l[i].Draw("same");
200  i++;
201  }
202 
203  TLatex tSubdet;
204  tSubdet.SetNDC();
205  tSubdet.SetTextAlign(21);
206  tSubdet.SetTextSize(0.027);
207  tSubdet.SetTextAngle(90);
208  for (unsigned int j=1;j<=6;j++ ){
209  auto thePart = static_cast<AlignmentPI::partitions>(j);
210  tSubdet.SetTextColor(kRed);
211  auto myPair = (j>1) ? AlignmentPI::calculatePosition(gPad,compare->GetBinLowEdge(boundaries[j-2])) : AlignmentPI::calculatePosition(gPad,compare->GetBinLowEdge(0));
212  float theX_ = myPair.first+0.025;
213  tSubdet.DrawLatex(theX_,0.20,Form("%s",(AlignmentPI::getStringFromPart(thePart)).c_str()));
214  }
215 
216  TLegend legend = TLegend(0.58,0.82,0.95,0.9);
217  legend.SetTextSize(0.03);
218  legend.SetHeader("Alignment comparison","C"); // option "C" allows to center the header
219  legend.AddEntry(compare.get(),("IOV:"+std::to_string(std::get<0>(lastiov))+"-"+std::to_string(std::get<0>(firstiov))).c_str(),"PL");
220  legend.Draw("same");
221 
222  TLatex t1;
223  t1.SetNDC();
224  t1.SetTextAlign(21);
225  t1.SetTextSize(0.05);
226  t1.DrawLatex(0.2, 0.93, Form("%s",s_coord.c_str()));
227  t1.SetTextColor(kBlue);
228  t1.DrawLatex(0.6, 0.93, Form("IOV %s - %s ",lastIOVsince.c_str(),firstIOVsince.c_str()));
229 
230  std::string fileName(m_imageFileName);
231  canvas.SaveAs(fileName.c_str());
232 
233  return true;
234  }
235  };
236 
237  typedef TrackerAlignmentCompare<AlignmentPI::t_x> TrackerAlignmentCompareX;
238  typedef TrackerAlignmentCompare<AlignmentPI::t_y> TrackerAlignmentCompareY;
239  typedef TrackerAlignmentCompare<AlignmentPI::t_z> TrackerAlignmentCompareZ;
240 
241  typedef TrackerAlignmentCompare<AlignmentPI::rot_alpha> TrackerAlignmentCompareAlpha;
242  typedef TrackerAlignmentCompare<AlignmentPI::rot_beta> TrackerAlignmentCompareBeta;
243  typedef TrackerAlignmentCompare<AlignmentPI::rot_gamma> TrackerAlignmentCompareGamma;
244 
245  //*******************************************//
246  // Summary canvas per subdetector
247  //******************************************//
248 
249  template<AlignmentPI::partitions q> class TrackerAlignmentSummary : public cond::payloadInspector::PlotImage<Alignments> {
250  public:
251  TrackerAlignmentSummary() : cond::payloadInspector::PlotImage<Alignments>( "Comparison of all coordinates between two geometries for "+getStringFromPart (q) ){
252  setSingleIov( false );
253  }
254 
255  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
256 
257  std::vector<std::tuple<cond::Time_t,cond::Hash> > sorted_iovs = iovs;
258 
259  // make absolute sure the IOVs are sortd by since
260  std::sort(begin(sorted_iovs), end(sorted_iovs), [](auto const &t1, auto const &t2) {
261  return std::get<0>(t1) < std::get<0>(t2);
262  });
263 
264  auto firstiov = sorted_iovs.front();
265  auto lastiov = sorted_iovs.back();
266 
267  std::shared_ptr<Alignments> last_payload = fetchPayload( std::get<1>(lastiov) );
268  std::shared_ptr<Alignments> first_payload = fetchPayload( std::get<1>(firstiov) );
269 
270  std::string lastIOVsince = std::to_string(std::get<0>(lastiov));
271  std::string firstIOVsince = std::to_string(std::get<0>(firstiov));
272 
273  std::vector<AlignTransform> ref_ali = first_payload->m_align;
274  std::vector<AlignTransform> target_ali = last_payload->m_align;
275 
276  if(ref_ali.size() != target_ali.size()){
277  edm::LogError("TrackerAlignment_PayloadInspector") << "the size of the reference alignment (" << ref_ali.size() << ") is different from the one of the target ("<< target_ali.size() << ")! You are probably trying to compare different underlying geometries. Exiting";
278  return false;
279  }
280 
281  // check that the geomtery is a tracker one
282  const char * path_toTopologyXML = (ref_ali.size()==AlignmentPI::phase0size) ? "Geometry/TrackerCommonData/data/trackerParameters.xml" : "Geometry/TrackerCommonData/data/PhaseI/trackerParameters.xml";
284 
285  for (const auto &ali : ref_ali){
286  auto mydetid = ali.rawId();
287  if(DetId(mydetid).det() != DetId::Tracker){
288  edm::LogWarning("TrackerAlignment_PayloadInspector") << "Encountered invalid Tracker DetId:" << DetId(mydetid).rawId() <<" ("<< DetId(mydetid).det()<<") is different from "<<DetId::Tracker<<" (is DoubleSide: "<< tTopo.tidIsDoubleSide(mydetid)<<"); subdetId "<< DetId(mydetid).subdetId() <<" - terminating ";
289  return false;
290  }
291  }
292 
293  TCanvas canvas("Alignment Comparison","Alignment Comparison",1800,1200);
294  canvas.Divide(3,2);
295 
296  std::unordered_map<AlignmentPI::coordinate,std::unique_ptr<TH1F> > diffs;
298 
299  for (const auto &coord : coords){
300 
301  auto s_coord = AlignmentPI::getStringFromCoordinate(coord);
302  std::string unit = (coord == AlignmentPI::t_x || coord == AlignmentPI::t_y || coord == AlignmentPI::t_z ) ? "[#mum]" : "[mrad]";
303 
304  diffs[coord] = std::make_unique<TH1F>(Form("hDiff_%s",s_coord.c_str()),Form(";#Delta%s %s;n. of modules",s_coord.c_str(),unit.c_str()),1000,-500.,500.);
305 
306  }
307 
308 
309  int loopedComponents(0);
310  for(unsigned int i=0;i<=ref_ali.size();i++){
311 
312  if(ref_ali[i].rawId() == target_ali[i].rawId()){
313 
314  loopedComponents++;
315  int subid = DetId(ref_ali[i].rawId()).subdetId();
316  auto thePart = static_cast<AlignmentPI::partitions>(subid);
317  if(thePart!=q) continue;
318 
319  CLHEP::HepRotation target_rot( target_ali[i].rotation() );
320  CLHEP::HepRotation ref_rot( ref_ali[i].rotation() );
321 
322  align::RotationType target_rotation( target_rot.xx(), target_rot.xy(), target_rot.xz(),
323  target_rot.yx(), target_rot.yy(), target_rot.yz(),
324  target_rot.zx(), target_rot.zy(), target_rot.zz() );
325 
326  align::RotationType ref_rotation( ref_rot.xx(), ref_rot.xy(), ref_rot.xz(),
327  ref_rot.yx(), ref_rot.yy(), ref_rot.yz(),
328  ref_rot.zx(), ref_rot.zy(), ref_rot.zz() );
329 
330  align::EulerAngles target_eulerAngles = align::toAngles(target_rotation);
331  align::EulerAngles ref_eulerAngles = align::toAngles(ref_rotation);
332 
333  for ( const auto &coord : coords){
334  switch(coord){
335  case AlignmentPI::t_x :
336  diffs[coord]->Fill((target_ali[i].translation().x()-ref_ali[i].translation().x())*AlignmentPI::cmToUm);
337  break;
338  case AlignmentPI::t_y:
339  diffs[coord]->Fill((target_ali[i].translation().y()-ref_ali[i].translation().y())*AlignmentPI::cmToUm);
340  break;
341  case AlignmentPI::t_z:
342  diffs[coord]->Fill((target_ali[i].translation().z()-ref_ali[i].translation().z())*AlignmentPI::cmToUm);
343  break;
345  diffs[coord]->Fill((target_eulerAngles[0]-ref_eulerAngles[0])*1000.);
346  break;
348  diffs[coord]->Fill((target_eulerAngles[1]-ref_eulerAngles[1])*1000.);
349  break;
350  case AlignmentPI::rot_gamma:
351  diffs[coord]->Fill((target_eulerAngles[2]-ref_eulerAngles[2])*1000.);
352  break;
353  default:
354  edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate "<< coord << std::endl;
355  break;
356  } // switch on the coordinate
357  }
358  } // check on the same detID
359  } // loop on the components
360 
361  int c_index=1;
362 
363  auto legend = std::unique_ptr<TLegend>(new TLegend(0.14,0.93,0.55,0.98));
364  legend->AddEntry(diffs[AlignmentPI::t_x].get(),("#DeltaIOV: "+std::to_string(std::get<0>(lastiov))+"-"+std::to_string(std::get<0>(firstiov))).c_str(),"L");
365  legend->SetTextSize(0.03);
366 
367  for (const auto &coord : coords){
368  canvas.cd(c_index)->SetLogy();
369  canvas.cd(c_index)->SetTopMargin(0.02);
370  canvas.cd(c_index)->SetBottomMargin(0.15);
371  canvas.cd(c_index)->SetLeftMargin(0.14);
372  canvas.cd(c_index)->SetRightMargin(0.04);
373  diffs[coord]->SetLineWidth(2);
374  AlignmentPI::makeNicePlotStyle(diffs[coord].get(),kBlack);
375 
376  //float x_max = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindLastBinAbove(0.));
377  //float x_min = diffs[coord]->GetXaxis()->GetBinCenter(diffs[coord]->FindFirstBinAbove(0.));
378  //float extremum = std::abs(x_max) > std::abs(x_min) ? std::abs(x_max) : std::abs(x_min);
379  //diffs[coord]->GetXaxis()->SetRangeUser(-extremum*2,extremum*2);
380 
381  int i_max = diffs[coord]->FindLastBinAbove(0.);
382  int i_min = diffs[coord]->FindFirstBinAbove(0.);
383  diffs[coord]->GetXaxis()->SetRange(std::max(1,i_min-10),std::min(i_max+10,diffs[coord]->GetNbinsX()));
384  diffs[coord]->Draw("HIST");
385  AlignmentPI::makeNiceStats(diffs[coord].get(),q,kBlack);
386 
387  legend->Draw("same");
388 
389  c_index++;
390  }
391 
392  std::string fileName(m_imageFileName);
393  canvas.SaveAs(fileName.c_str());
394 
395  return true;
396  }
397  };
398 
399  typedef TrackerAlignmentSummary<AlignmentPI::BPix> TrackerAlignmentSummaryBPix;
400  typedef TrackerAlignmentSummary<AlignmentPI::FPix> TrackerAlignmentSummaryFPix;
401  typedef TrackerAlignmentSummary<AlignmentPI::TIB> TrackerAlignmentSummaryTIB;
402 
403  typedef TrackerAlignmentSummary<AlignmentPI::TID> TrackerAlignmentSummaryTID;
404  typedef TrackerAlignmentSummary<AlignmentPI::TOB> TrackerAlignmentSummaryTOB;
405  typedef TrackerAlignmentSummary<AlignmentPI::TEC> TrackerAlignmentSummaryTEC;
406 
407 
408  //*******************************************//
409  // History of the position of the BPix Barycenter
410  //******************************************//
411 
412  template<AlignmentPI::coordinate coord> class BPixBarycenterHistory : public cond::payloadInspector::HistoryPlot<Alignments,float> {
413  public:
414  BPixBarycenterHistory(): cond::payloadInspector::HistoryPlot<Alignments,float>(" Barrel Pixel "+AlignmentPI::getStringFromCoordinate(coord)+" positions vs time", AlignmentPI::getStringFromCoordinate(coord)+" position [cm]"){}
415  ~BPixBarycenterHistory() override = default;
416 
417  float getFromPayload( Alignments& payload ) override {
418 
419  std::vector<AlignTransform> alignments = payload.m_align;
420 
421  float barycenter=0.;
422  float nmodules(0.);
423  for(const auto& ali : alignments ){
424 
425  if(DetId(ali.rawId()).det() != DetId::Tracker){
426  edm::LogWarning("TrackerAlignment_PayloadInspector") << "Encountered invalid Tracker DetId:" << ali.rawId() <<" "<<DetId(ali.rawId()).det()<<" is different from "<<DetId::Tracker<<" - terminating ";
427  return false;
428  }
429 
430  int subid = DetId(ali.rawId()).subdetId();
431  if(subid!=PixelSubdetector::PixelBarrel) continue;
432 
433  nmodules++;
434  switch(coord){
435  case AlignmentPI::t_x :
436  barycenter+=(ali.translation().x());
437  break;
438  case AlignmentPI::t_y:
439  barycenter+=(ali.translation().y());
440  break;
441  case AlignmentPI::t_z:
442  barycenter+=(ali.translation().z());
443  break;
444  default:
445  edm::LogError("TrackerAlignment_PayloadInspector") << "Unrecognized coordinate "<< coord << std::endl;
446  break;
447  } // switch on the coordinate (only X,Y,Z are interesting)
448  } // ends loop on the alignments
449 
450  edm::LogInfo("TrackerAlignment_PayloadInspector")<<"barycenter ("<<barycenter<<")/n. modules ("<< nmodules << ") = "<< barycenter/nmodules << std::endl;
451 
452  // take the mean
453  barycenter/=nmodules;
454 
455  // applied GPR correction to move barycenter to global CMS coordinates
456  barycenter+=hardcodeGPR.at(coord);
457 
458  return barycenter;
459 
460  } // payload
461  };
462 
463  typedef BPixBarycenterHistory<AlignmentPI::t_x> X_BPixBarycenterHistory;
464  typedef BPixBarycenterHistory<AlignmentPI::t_y> Y_BPixBarycenterHistory;
465  typedef BPixBarycenterHistory<AlignmentPI::t_z> Z_BPixBarycenterHistory;
466 
467 
468 }
469 
471  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareX);
472  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareY);
473  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareZ);
474  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareAlpha);
475  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareBeta);
476  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentCompareGamma);
477  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryBPix);
478  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryFPix);
479  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTIB);
480  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTID);
481  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTOB);
482  PAYLOAD_INSPECTOR_CLASS(TrackerAlignmentSummaryTEC);
483  PAYLOAD_INSPECTOR_CLASS(X_BPixBarycenterHistory);
484  PAYLOAD_INSPECTOR_CLASS(Y_BPixBarycenterHistory);
485  PAYLOAD_INSPECTOR_CLASS(Z_BPixBarycenterHistory);
486 }
void makeNicePlotStyle(TH1 *hist, int color)
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
std::string getStringFromPart(AlignmentPI::partitions i)
std::vector< AlignTransform > m_align
Definition: Alignments.h:19
static const float cmToUm
void makeNiceStats(TH1F *hist, AlignmentPI::partitions part, int color)
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
static const unsigned int phase0size
std::pair< double, double > calculatePosition(TVirtualPad *myPad, int boundary)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
virtual bool fill(const std::vector< std::tuple< cond::Time_t, cond::Hash > > &iovs)=0
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define end
Definition: vmac.h:39
T min(T a, T b)
Definition: MathUtil.h:58
EulerAngles toAngles(const RotationType &)
Convert rotation matrix to angles about x-, y-, z-axes (frame rotation).
Definition: Utilities.cc:9
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Definition: DetId.h:18
AlgebraicVector EulerAngles
Definition: Definitions.h:36
std::string getStringFromCoordinate(AlignmentPI::coordinate coord)
bool tidIsDoubleSide(const DetId &id) const
Definition: plugin.cc:24
#define begin
Definition: vmac.h:32
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
def canvas(sub, attr)
Definition: svgfig.py:482
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39