CMS 3D CMS Logo

SiPixelBarycenterHelper.h
Go to the documentation of this file.
1 #ifndef DQM_SiPixelPhase1Summary_SiPixelBarycenterHelper_h
2 #define DQM_SiPixelPhase1Summary_SiPixelBarycenterHelper_h
3 
10 
11 // Helper mainly based on https://github.com/cms-sw/cmssw/blob/master/CondCore/AlignmentPlugins/interface/AlignmentPayloadInspectorHelper.h
12 
13 namespace DQMBarycenter {
14 
15  enum coordinate {
16  t_x = 1,
17  t_y = 2,
18  t_z = 3,
19  rot_alpha = 4,
20  rot_beta = 5,
21  rot_gamma = 6,
22  };
23 
24  enum partitions {
25  BPix = 1,
26  FPix_zp = 2,
27  FPix_zm = 3,
28  BPix_xp = 4,
29  BPix_xm = 5,
34  };
35 
36  enum class PARTITION {
37  BPIX, // 0 Barrel Pixel
38  FPIX_zp, // 1 Forward Pixel Z-Plus
39  FPIX_zm, // 2 Forward Pixel Z-Minus
40  BPIX_xp, // 3 Barrel Pixel X-Plus
41  BPIX_xm, // 4 Barrel Pixel X-Minus
42  FPIX_zp_xp, // 5 Forward Pixel Z-Plus X-Plus
43  FPIX_zm_xp, // 6 Forward Pixel Z-Minus X-Plus
44  FPIX_zp_xm, // 7 Forward Pixel Z-Plus X-Minus
45  FPIX_zm_xm, // 8 Forward Pixel Z-Minus X-Minus
47  };
48 
49  extern const PARTITION PARTITIONS[(int)PARTITION::LAST + 1];
50  const PARTITION PARTITIONS[] = {
60  };
61 
63  std::map<DQMBarycenter::PARTITION, double> Xbarycenters;
64  std::map<DQMBarycenter::PARTITION, double> Ybarycenters;
65  std::map<DQMBarycenter::PARTITION, double> Zbarycenters;
66  std::map<DQMBarycenter::PARTITION, double> nmodules;
67 
68  public:
69  inline void init();
70 
71  /*--------------------------------------------------------------------*/
72  inline const std::array<double, 9> getX()
73  /*--------------------------------------------------------------------*/
74  {
75  return {{Xbarycenters[PARTITION::BPIX],
84  };
85 
86  /*--------------------------------------------------------------------*/
87  inline const std::array<double, 9> getY()
88  /*--------------------------------------------------------------------*/
89  {
90  return {{Ybarycenters[PARTITION::BPIX],
99  };
100 
101  /*--------------------------------------------------------------------*/
102  inline const std::array<double, 9> getZ()
103  /*--------------------------------------------------------------------*/
104  {
105  return {{Zbarycenters[PARTITION::BPIX],
114  };
115  inline virtual ~TkAlBarycenters() {}
116 
117  /*--------------------------------------------------------------------*/
118  inline void computeBarycenters(const std::vector<AlignTransform>& input,
119  const TrackerTopology& tTopo,
120  const std::map<DQMBarycenter::coordinate, float>& GPR)
121  /*--------------------------------------------------------------------*/
122  {
123  for (const auto& ali : input) {
124  if (DetId(ali.rawId()).det() != DetId::Tracker) {
125  edm::LogWarning("SiPixelBarycenters::computeBarycenters")
126  << "Encountered invalid Tracker DetId:" << ali.rawId() << " " << DetId(ali.rawId()).det()
127  << " is different from " << DetId::Tracker << " - terminating ";
128  assert(DetId(ali.rawId()).det() != DetId::Tracker);
129  }
130 
131  const auto& tns = align::TrackerNameSpace(&tTopo);
132  int subid = DetId(ali.rawId()).subdetId();
133  if (subid == PixelSubdetector::PixelBarrel || subid == PixelSubdetector::PixelEndcap) { // use only pixel
134  switch (subid) { // Separate BPIX, FPIX_zp and FPIX_zm
136  Xbarycenters[PARTITION::BPIX] += (ali.translation().x());
137  Ybarycenters[PARTITION::BPIX] += (ali.translation().y());
138  Zbarycenters[PARTITION::BPIX] += (ali.translation().z());
140  break;
142  // minus side
143  if (tns.tpe().endcapNumber(ali.rawId()) == 1) {
144  Xbarycenters[PARTITION::FPIX_zm] += (ali.translation().x());
145  Ybarycenters[PARTITION::FPIX_zm] += (ali.translation().y());
146  Zbarycenters[PARTITION::FPIX_zm] += (ali.translation().z());
148  } // plus side
149  else {
150  Xbarycenters[PARTITION::FPIX_zp] += (ali.translation().x());
151  Ybarycenters[PARTITION::FPIX_zp] += (ali.translation().y());
152  Zbarycenters[PARTITION::FPIX_zp] += (ali.translation().z());
154  }
155  break;
156  default:
157  edm::LogError("PixelDQM") << "Unrecognized partition for barycenter computation " << subid << std::endl;
158  break;
159  }
160 
161  switch (subid) { // Separate following the PCL HLS
163  if ((PixelBarrelName(DetId(ali.rawId()), true).shell() == PixelBarrelName::mO) ||
164  (PixelBarrelName(DetId(ali.rawId()), true).shell() == PixelBarrelName::pO)) { // BPIX x-
165  Xbarycenters[PARTITION::BPIX_xm] += (ali.translation().x());
166  Ybarycenters[PARTITION::BPIX_xm] += (ali.translation().y());
167  Zbarycenters[PARTITION::BPIX_xm] += (ali.translation().z());
169  } else { // BPIX x+
170  Xbarycenters[PARTITION::BPIX_xp] += (ali.translation().x());
171  Ybarycenters[PARTITION::BPIX_xp] += (ali.translation().y());
172  Zbarycenters[PARTITION::BPIX_xp] += (ali.translation().z());
174  }
175  break;
177  if (PixelEndcapName(DetId(ali.rawId()), true).halfCylinder() == PixelEndcapName::mO) { //FPIX z- x-
178  Xbarycenters[PARTITION::FPIX_zm_xm] += (ali.translation().x());
179  Ybarycenters[PARTITION::FPIX_zm_xm] += (ali.translation().y());
180  Zbarycenters[PARTITION::FPIX_zm_xm] += (ali.translation().z());
182  } else if (PixelEndcapName(DetId(ali.rawId()), true).halfCylinder() ==
183  PixelEndcapName::mI) { //FPIX z- x+
184  Xbarycenters[PARTITION::FPIX_zm_xp] += (ali.translation().x());
185  Ybarycenters[PARTITION::FPIX_zm_xp] += (ali.translation().y());
186  Zbarycenters[PARTITION::FPIX_zm_xp] += (ali.translation().z());
188  } else if (PixelEndcapName(DetId(ali.rawId()), true).halfCylinder() ==
189  PixelEndcapName::pO) { //FPIX z+ x-
190  Xbarycenters[PARTITION::FPIX_zp_xm] += (ali.translation().x());
191  Ybarycenters[PARTITION::FPIX_zp_xm] += (ali.translation().y());
192  Zbarycenters[PARTITION::FPIX_zp_xm] += (ali.translation().z());
194  } else if (PixelEndcapName(DetId(ali.rawId()), true).halfCylinder() ==
195  PixelEndcapName::pI) { //FPIX z+ x+
196  Xbarycenters[PARTITION::FPIX_zp_xp] += (ali.translation().x());
197  Ybarycenters[PARTITION::FPIX_zp_xp] += (ali.translation().y());
198  Zbarycenters[PARTITION::FPIX_zp_xp] += (ali.translation().z());
200  } else {
201  edm::LogError("PixelDQM") << "Unrecognized partition for barycenter computation " << subid << std::endl;
202  }
203  break;
204  default:
205  edm::LogError("PixelDQM") << "Unrecognized partition for barycenter computation " << subid << std::endl;
206  break;
207  }
208  }
209  }
210 
211  for (const auto& p : PARTITIONS) {
212  Xbarycenters[p] /= nmodules[p];
213  Ybarycenters[p] /= nmodules[p];
214  Zbarycenters[p] /= nmodules[p];
215 
216  Xbarycenters[p] += GPR.at(DQMBarycenter::t_x);
217  Ybarycenters[p] += GPR.at(DQMBarycenter::t_y);
218  Zbarycenters[p] += GPR.at(DQMBarycenter::t_z);
219  }
220  }
221  };
222 } // namespace DQMBarycenter
223 
224 #endif
DQMBarycenter::TkAlBarycenters::getX
const std::array< double, 9 > getX()
Definition: SiPixelBarycenterHelper.h:72
PixelEndcapName::mO
Definition: PixelEndcapName.h:18
DQMBarycenter::PARTITION::BPIX_xp
PixelBarrelName::pO
Definition: PixelBarrelName.h:18
input
static const std::string input
Definition: EdmProvDump.cc:48
PixelBarrelName.h
MessageLogger.h
PixelSubdetector::PixelEndcap
Definition: PixelSubdetector.h:11
PixelSubdetector::PixelBarrel
Definition: PixelSubdetector.h:11
DQMBarycenter::PARTITION
PARTITION
Definition: SiPixelBarycenterHelper.h:36
DQMBarycenter::PARTITION::BPIX_xm
TrackerTopology
Definition: TrackerTopology.h:16
DQMBarycenter::TkAlBarycenters::Xbarycenters
std::map< DQMBarycenter::PARTITION, double > Xbarycenters
Definition: SiPixelBarycenterHelper.h:63
DQMBarycenter::PARTITION::FPIX_zp_xp
DQMBarycenter::FPix_zp_xm
Definition: SiPixelBarycenterHelper.h:32
PixelBarrelName
Definition: PixelBarrelName.h:16
Alignments.h
cms::cuda::assert
assert(be >=bs)
DQMBarycenter::FPix_zp_xp
Definition: SiPixelBarycenterHelper.h:30
DQMBarycenter::t_x
Definition: SiPixelBarycenterHelper.h:16
DQMBarycenter::FPix_zm_xm
Definition: SiPixelBarycenterHelper.h:33
PixelEndcapName::halfCylinder
HalfCylinder halfCylinder() const
Definition: PixelEndcapName.h:42
DQMBarycenter::TkAlBarycenters::getY
const std::array< double, 9 > getY()
Definition: SiPixelBarycenterHelper.h:87
PixelBarrelName::mO
Definition: PixelBarrelName.h:18
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
align::TrackerNameSpace
Definition: TrackerNameSpace.h:10
DQMBarycenter::TkAlBarycenters::nmodules
std::map< DQMBarycenter::PARTITION, double > nmodules
Definition: SiPixelBarycenterHelper.h:66
DQMBarycenter::BPix_xm
Definition: SiPixelBarycenterHelper.h:29
DQMBarycenter::rot_beta
Definition: SiPixelBarycenterHelper.h:20
DQMBarycenter::PARTITIONS
const PARTITION PARTITIONS[(int) PARTITION::LAST+1]
Definition: SiPixelBarycenterHelper.h:50
TrackerNameSpace.h
DetId
Definition: DetId.h:17
TrackerTopology.h
PixelEndcapName
Definition: PixelEndcapName.h:16
DQMBarycenter::BPix_xp
Definition: SiPixelBarycenterHelper.h:28
DQMBarycenter::partitions
partitions
Definition: SiPixelBarycenterHelper.h:24
DQMBarycenter::PARTITION::FPIX_zm
DQMBarycenter::PARTITION::BPIX
DQMBarycenter::t_y
Definition: SiPixelBarycenterHelper.h:17
PixelEndcapName::mI
Definition: PixelEndcapName.h:18
DQMBarycenter::BPix
Definition: SiPixelBarycenterHelper.h:25
DQMBarycenter::TkAlBarycenters::Ybarycenters
std::map< DQMBarycenter::PARTITION, double > Ybarycenters
Definition: SiPixelBarycenterHelper.h:64
DetId::Tracker
Definition: DetId.h:25
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
DQMBarycenter::PARTITION::FPIX_zp_xm
PixelBarrelName::shell
Shell shell() const
Definition: PixelBarrelName.h:40
DQMBarycenter::PARTITION::LAST
DQMBarycenter::TkAlBarycenters::getZ
const std::array< double, 9 > getZ()
Definition: SiPixelBarycenterHelper.h:102
DQMBarycenter::FPix_zm_xp
Definition: SiPixelBarycenterHelper.h:31
createfilelist.int
int
Definition: createfilelist.py:10
DQMBarycenter::FPix_zp
Definition: SiPixelBarycenterHelper.h:26
DQMBarycenter::PARTITION::FPIX_zm_xm
DQMBarycenter::PARTITION::FPIX_zm_xp
DQMBarycenter::t_z
Definition: SiPixelBarycenterHelper.h:18
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
DQMBarycenter::coordinate
coordinate
Definition: SiPixelBarycenterHelper.h:15
DQMBarycenter::TkAlBarycenters
Definition: SiPixelBarycenterHelper.h:62
PixelEndcapName::pO
Definition: PixelEndcapName.h:18
DQMBarycenter
Definition: SiPixelBarycenterHelper.h:13
DQMBarycenter::TkAlBarycenters::computeBarycenters
void computeBarycenters(const std::vector< AlignTransform > &input, const TrackerTopology &tTopo, const std::map< DQMBarycenter::coordinate, float > &GPR)
Definition: SiPixelBarycenterHelper.h:118
DQMBarycenter::TkAlBarycenters::~TkAlBarycenters
virtual ~TkAlBarycenters()
Definition: SiPixelBarycenterHelper.h:115
DQMBarycenter::FPix_zm
Definition: SiPixelBarycenterHelper.h:27
DQMBarycenter::TkAlBarycenters::Zbarycenters
std::map< DQMBarycenter::PARTITION, double > Zbarycenters
Definition: SiPixelBarycenterHelper.h:65
DQMBarycenter::TkAlBarycenters::init
void init()
DQMBarycenter::rot_alpha
Definition: SiPixelBarycenterHelper.h:19
PixelEndcapName::pI
Definition: PixelEndcapName.h:18
DQMBarycenter::PARTITION::FPIX_zp
DQMBarycenter::rot_gamma
Definition: SiPixelBarycenterHelper.h:21
PixelEndcapName.h