CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
HGCalCellUV Class Reference

#include <HGCalCellUV.h>

Public Member Functions

std::pair< int32_t, int32_t > cellUVFromXY1 (double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const
 
std::pair< int32_t, int32_t > cellUVFromXY1 (double xloc, double yloc, int32_t placement, int32_t type, int32_t partial, bool extend, bool debug) const
 
std::pair< int32_t, int32_t > cellUVFromXY2 (double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const
 
std::pair< int32_t, int32_t > cellUVFromXY2 (double xloc, double yloc, int32_t placement, int32_t type, int32_t partial, bool extend, bool debug) const
 
std::pair< int32_t, int32_t > cellUVFromXY3 (double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const
 
std::pair< int32_t, int32_t > cellUVFromXY4 (double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug)
 
 HGCalCellUV (double waferSize, double separation, int32_t nFine, int32_t nCoarse)
 

Private Member Functions

std::pair< int32_t, int32_t > cellUVFromXY4 (double xloc, double yloc, int ncell, double cellX, double cellY, double cellXTotal, double cellYTotal, std::map< std::pair< int, int >, std::pair< double, double > > &cellPos, bool extend, bool debug)
 

Private Attributes

std::map< std::pair< int32_t, int32_t >, std::pair< double, double > > cellPosCoarse_ [HGCalCell::cellPlacementTotal]
 
std::map< std::pair< int32_t, int32_t >, std::pair< double, double > > cellPosFine_ [HGCalCell::cellPlacementTotal]
 
double cellX_ [2]
 
double cellXTotal_ [2]
 
double cellY_ [2]
 
double cellYTotal_ [2]
 
std::unique_ptr< HGCalCellhgcalcell_
 
int32_t ncell_ [2]
 
const double waferSize_
 

Static Private Attributes

static constexpr double cos60_ = 0.5
 
static constexpr double sin60_ = 0.5 * sqrt3_
 
static constexpr double sqrt3_ = 1.732050807568877
 

Detailed Description

Definition at line 11 of file HGCalCellUV.h.

Constructor & Destructor Documentation

◆ HGCalCellUV()

HGCalCellUV::HGCalCellUV ( double  waferSize,
double  separation,
int32_t  nFine,
int32_t  nCoarse 
)

Definition at line 10 of file HGCalCellUV.cc.

References cms::cuda::assert(), HGCalCell::cellPlacementExtra, HGCalCell::cellPlacementTotal, cellPosCoarse_, cellPosFine_, cellX_, cellXTotal_, cellY_, cellYTotal_, hgcalcell_, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::iv, dqmdumpme::k, ncell_, sqrt3_, and findQualityFiles::v.

10  : waferSize_(waferSize) {
11  hgcalcell_ = std::make_unique<HGCalCell>(waferSize, nFine, nCoarse);
12  assert(nFine > 0 && nCoarse > 0);
13  ncell_[0] = nFine;
14  ncell_[1] = nCoarse;
15  for (int k = 0; k < 2; ++k) {
16  cellX_[k] = waferSize / (3 * ncell_[k]);
17  cellY_[k] = 0.5 * sqrt3_ * cellX_[k];
18  cellXTotal_[k] = (waferSize + separation) / (3 * ncell_[k]);
19  cellYTotal_[k] = 0.5 * sqrt3_ * cellXTotal_[k];
20  }
21  // Fill up the placement vectors
22  for (int placement = 0; placement < HGCalCell::cellPlacementTotal; ++placement) {
23  // Fine cells
24  for (int iu = 0; iu < 2 * ncell_[0]; ++iu) {
25  for (int iv = 0; iv < 2 * ncell_[0]; ++iv) {
26  int u = (placement < HGCalCell::cellPlacementExtra) ? iv : iu;
27  int v = (placement < HGCalCell::cellPlacementExtra) ? iu : iv;
28  if (((v - u) < ncell_[0]) && (u - v) <= ncell_[0]) {
29  cellPosFine_[placement][std::pair<int, int>(u, v)] = hgcalcell_->cellUV2XY1(u, v, placement, 0);
30  }
31  }
32  }
33  // Coarse cells
34  for (int iu = 0; iu < 2 * ncell_[1]; ++iu) {
35  for (int iv = 0; iv < 2 * ncell_[1]; ++iv) {
36  int u = (placement < HGCalCell::cellPlacementExtra) ? iv : iu;
37  int v = (placement < HGCalCell::cellPlacementExtra) ? iu : iv;
38  if (((v - u) < ncell_[1]) && (u - v) <= ncell_[1]) {
39  cellPosCoarse_[placement][std::pair<int, int>(u, v)] = hgcalcell_->cellUV2XY1(u, v, placement, 1);
40  }
41  }
42  }
43  }
44 }
static constexpr double sqrt3_
Definition: HGCalCellUV.h:51
std::map< std::pair< int32_t, int32_t >, std::pair< double, double > > cellPosFine_[HGCalCell::cellPlacementTotal]
Definition: HGCalCellUV.h:59
static constexpr int32_t cellPlacementTotal
Definition: HGCalCell.h:26
double cellXTotal_[2]
Definition: HGCalCellUV.h:57
double cellYTotal_[2]
Definition: HGCalCellUV.h:57
assert(be >=bs)
double cellY_[2]
Definition: HGCalCellUV.h:57
double cellX_[2]
Definition: HGCalCellUV.h:57
const double waferSize_
Definition: HGCalCellUV.h:55
std::unique_ptr< HGCalCell > hgcalcell_
Definition: HGCalCellUV.h:58
int32_t ncell_[2]
Definition: HGCalCellUV.h:56
std::map< std::pair< int32_t, int32_t >, std::pair< double, double > > cellPosCoarse_[HGCalCell::cellPlacementTotal]
Definition: HGCalCellUV.h:59
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

Member Function Documentation

◆ cellUVFromXY1() [1/2]

std::pair< int32_t, int32_t > HGCalCellUV::cellUVFromXY1 ( double  xloc,
double  yloc,
int32_t  placement,
int32_t  type,
bool  extend,
bool  debug 
) const

Definition at line 46 of file HGCalCellUV.cc.

References HGCalCell::cellPlacementExtra, cellY_, cellYTotal_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), cos60_, debug, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::iv, ncell_, makeMuonMisalignmentScenario::rot, sin60_, findQualityFiles::v, w(), x, and y.

Referenced by HGCalMouseBiteTester::analyze(), cellUVFromXY1(), and cellUVFromXY2().

47  {
48  if (type != 0)
49  type = 1;
50  //--- Reverse transform to placement=0, if placement index ≠ 6
51  double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
52  double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
53  int rot = placement % HGCalCell::cellPlacementExtra;
54  static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
55  static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
56  double x = xloc1 * fcos[rot] - yloc * fsin[rot];
57  double y = xloc1 * fsin[rot] + yloc * fcos[rot];
58 
59  //--- Calculate coordinates in u,v,w system
60  double u = x * sin60_ + y * cos60_;
61  double v = -x * sin60_ + y * cos60_;
62  double w = y;
63 
64  //--- Rounding in u, v, w coordinates
65  int iu = std::floor(u / cellY) + 3 * (ncell_[type]) - 1;
66  int iv = std::floor(v / cellY) + 3 * (ncell_[type]);
67  int iw = std::floor(w / cellY) + 1;
68 
69  int isv = (iu + iw) / 3;
70  int isu = (iv + iw) / 3;
71 
72  //--- Taking care of extending cells
73  if ((iu + iw) < 0) {
74  isu = (iv + iw + 1) / 3;
75  isv = 0;
76  } else if (isv - isu > ncell_[type] - 1) {
77  isu = (iv + iw + 1) / 3;
78  isv = (iu + iw - 1) / 3;
79  } else if (isu > 2 * ncell_[type] - 1) {
80  isu = 2 * ncell_[type] - 1;
81  isv = (iu + iw - 1) / 3;
82  }
83  if (debug)
84  edm::LogVerbatim("HGCalGeom") << "cellUVFromXY1: Input " << xloc << ":" << yloc << ":" << extend << " Output "
85  << isu << ":" << isv;
86  return std::make_pair(isu, isv);
87 }
Log< level::Info, true > LogVerbatim
T w() const
double cellYTotal_[2]
Definition: HGCalCellUV.h:57
double cellY_[2]
Definition: HGCalCellUV.h:57
static constexpr double cos60_
Definition: HGCalCellUV.h:53
#define debug
Definition: HDRShower.cc:19
int32_t ncell_[2]
Definition: HGCalCellUV.h:56
static constexpr double sin60_
Definition: HGCalCellUV.h:52
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

◆ cellUVFromXY1() [2/2]

std::pair< int32_t, int32_t > HGCalCellUV::cellUVFromXY1 ( double  xloc,
double  yloc,
int32_t  placement,
int32_t  type,
int32_t  partial,
bool  extend,
bool  debug 
) const

Definition at line 281 of file HGCalCellUV.cc.

References HGCalCell::cellPlacementExtra, cellUVFromXY1(), cellX_, cellXTotal_, cellY_, cellYTotal_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), cos60_, debug, HGCalTypes::edgeWaferHDBottom, HGCalTypes::edgeWaferLDTop, SiStripPI::max, ncell_, makeMuonMisalignmentScenario::rot, sin60_, sqrt3_, findQualityFiles::v, HGCalTypes::WaferHDBottom, and HGCalTypes::WaferLDTop.

288  {
289  if (type != 0)
290  type = 1;
291  double cellX = (extend) ? cellXTotal_[type] : cellX_[type];
292  double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
293  std::pair<int, int> uv = HGCalCellUV::cellUVFromXY1(xloc, yloc, placement, type, extend, debug);
294  int u = uv.first;
295  int v = uv.second;
296  if (partial == HGCalTypes::WaferLDTop) {
298  double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
299  int rot = placement % HGCalCell::cellPlacementExtra;
300  static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
301  static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
302  double xprime = -1 * (xloc1 * fcos[rot] - yloc * fsin[rot]);
303  double yprime = xloc1 * fsin[rot] + yloc * fcos[rot];
304  double xcell = -1 * (1.5 * (v - u) + 0.5) * cellX;
305  double ycell = (v + u - 2 * ncell_[type] + 1) * cellY;
306  if ((yprime - sqrt3_ * xprime) > (ycell - sqrt3_ * xcell)) {
307  --u;
308  if ((v - u) >= ncell_[1])
309  --v;
310  } else {
311  --u;
312  --v;
313  v = std::max(v, 0);
314  }
315  }
316  } else if (partial == HGCalTypes::WaferHDBottom) {
319  double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
320  int rot = placement % HGCalCell::cellPlacementExtra;
321  static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
322  static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
323  double xprime = -1 * (xloc1 * fcos[rot] - yloc * fsin[rot]);
324  double yprime = xloc1 * fsin[rot] + yloc * fcos[rot];
325  double xcell = -1 * (1.5 * (v - u) + 0.5) * cellX;
326  double ycell = (v + u - 2 * ncell_[type] + 1) * cellY;
327  if ((yprime - sqrt3_ * xprime) > (ycell - sqrt3_ * xcell)) {
328  ++u;
329  ++v;
330  } else {
331  ++u;
332  }
333  }
334  }
335  if (debug)
336  edm::LogVerbatim("HGCalGeom") << "cellUVFromXY5: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
337  << ":" << v;
338  return std::make_pair(u, v);
339 }
static constexpr std::array< int, 3 > edgeWaferHDBottom
Definition: HGCalTypes.h:88
Log< level::Info, true > LogVerbatim
static constexpr std::array< int, 3 > edgeWaferLDTop
Definition: HGCalTypes.h:81
std::pair< int32_t, int32_t > cellUVFromXY1(double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const
Definition: HGCalCellUV.cc:46
static constexpr double sqrt3_
Definition: HGCalCellUV.h:51
double cellXTotal_[2]
Definition: HGCalCellUV.h:57
double cellYTotal_[2]
Definition: HGCalCellUV.h:57
double cellY_[2]
Definition: HGCalCellUV.h:57
double cellX_[2]
Definition: HGCalCellUV.h:57
static constexpr int32_t WaferHDBottom
Definition: HGCalTypes.h:52
static constexpr double cos60_
Definition: HGCalCellUV.h:53
#define debug
Definition: HDRShower.cc:19
int32_t ncell_[2]
Definition: HGCalCellUV.h:56
static constexpr double sin60_
Definition: HGCalCellUV.h:52
static constexpr int32_t WaferLDTop
Definition: HGCalTypes.h:45
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

◆ cellUVFromXY2() [1/2]

std::pair< int32_t, int32_t > HGCalCellUV::cellUVFromXY2 ( double  xloc,
double  yloc,
int32_t  placement,
int32_t  type,
bool  extend,
bool  debug 
) const

Definition at line 89 of file HGCalCellUV.cc.

References alignmentValidation::c1, reco::ceil(), HGCalCell::cellPlacementExtra, cellY_, cellYTotal_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), cos60_, debug, ncell_, alignCSCRings::r, makeMuonMisalignmentScenario::rot, sin60_, sqrt3_, testProducerWithPsetDescEmpty_cfi::u1, MetAnalyzer::u2, testProducerWithPsetDescEmpty_cfi::u3, interactiveExample::ui, findQualityFiles::v, x, and y.

Referenced by HGCalMouseBiteTester::analyze().

90  {
91  //--- Using multiple inequalities to find (u, v)
92  //--- Reverse transform to placement=0, if placement index ≠ 7
93  if (type != 0)
94  type = 1;
95  double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
96  double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -1 * xloc;
97  int rot = placement % HGCalCell::cellPlacementExtra;
98  static constexpr std::array<double, 6> fcos = {{cos60_, 1.0, cos60_, -cos60_, -1.0, -cos60_}};
99  static constexpr std::array<double, 6> fsin = {{-sin60_, 0.0, sin60_, sin60_, 0.0, -sin60_}};
100  double x = xloc1 * fcos[rot] - yloc * fsin[rot];
101  double y = xloc1 * fsin[rot] + yloc * fcos[rot];
102 
103  int32_t u(-100), v(-100);
104  int ncell = ncell_[type];
105  double r = cellY;
106  double l1 = (y / r) + ncell - 1.0;
107  int l2 = std::floor((0.5 * y + 0.5 * x / sqrt3_) / r + ncell - 4.0 / 3.0);
108  int l3 = std::floor((x / sqrt3_) / r + ncell - 4.0 / 3.0);
109  double l4 = (y + sqrt3_ * x) / (2 * r) + 2 * ncell - 2;
110  double l5 = (y - sqrt3_ * x) / (2 * r) - ncell;
111  double u1 = (y / r) + ncell + 1.0;
112  int u2 = std::ceil((0.5 * y + 0.5 * x / sqrt3_) / r + ncell + 2.0 / 3.0);
113  int u3 = std::ceil((x / sqrt3_) / r + ncell);
114  double u4 = l4 + 2;
115  double u5 = l5 + 2;
116 
117  for (int ui = l2 + 1; ui < u2; ui++) {
118  for (int vi = l3 + 1; vi < u3; vi++) {
119  int c1 = 2 * ui - vi;
120  int c4 = ui + vi;
121  int c5 = ui - 2 * vi;
122  if ((c1 < u1) && (c1 > l1) && (c4 < u4) && (c4 > l4) && (c5 < u5) && (c5 > l5)) {
123  u = ui;
124  v = vi;
125  }
126  }
127  }
128 
129  //--- Taking care of extending cells
130  if (v == -1) {
131  if (y < (2 * u - v - ncell) * r) {
132  ++v;
133  } else {
134  ++u;
135  ++v;
136  }
137  }
138  if (v - u == ncell) {
139  if ((y + sqrt3_ * x) < ((u + v - 2 * ncell + 1) * 2 * r)) {
140  --v;
141  } else {
142  ++u;
143  }
144  }
145  if (u == 2 * ncell) {
146  if ((y - sqrt3_ * x) < ((u - 2 * v + ncell - 1) * 2 * r)) {
147  --u;
148  } else {
149  --u;
150  --v;
151  }
152  }
153  if (debug)
154  edm::LogVerbatim("HGCalGeom") << "cellUVFromXY2: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
155  << ":" << v;
156  return std::make_pair(u, v);
157 }
Log< level::Info, true > LogVerbatim
constexpr int32_t ceil(float num)
static constexpr double sqrt3_
Definition: HGCalCellUV.h:51
double cellYTotal_[2]
Definition: HGCalCellUV.h:57
double cellY_[2]
Definition: HGCalCellUV.h:57
static constexpr double cos60_
Definition: HGCalCellUV.h:53
#define debug
Definition: HDRShower.cc:19
int32_t ncell_[2]
Definition: HGCalCellUV.h:56
static constexpr double sin60_
Definition: HGCalCellUV.h:52
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

◆ cellUVFromXY2() [2/2]

std::pair< int32_t, int32_t > HGCalCellUV::cellUVFromXY2 ( double  xloc,
double  yloc,
int32_t  placement,
int32_t  type,
int32_t  partial,
bool  extend,
bool  debug 
) const

Definition at line 341 of file HGCalCellUV.cc.

References HGCalCell::cellPlacementExtra, cellUVFromXY1(), debug, HGCalTypes::edgeWaferHDBottom, HGCalTypes::edgeWaferLDTop, hgcalcell_, HGCalWaferMask::maskCut(), SiStripPI::max, ncell_, findQualityFiles::v, HGCalTypes::WaferFull, HGCalTypes::WaferHDBottom, HGCalTypes::WaferHDTop, HGCalTypes::WaferLDThree, HGCalTypes::WaferLDTop, and waferSize_.

348  {
349  if (type != 0)
350  type = 1;
351  std::pair<int, int> uv = HGCalCellUV::cellUVFromXY1(xloc, yloc, placement, type, extend, debug);
352  int u = uv.first;
353  int v = uv.second;
354  if ((partial != HGCalTypes::WaferFull) && (type == 1)) {
355  if (u == 1 && v == 8) {
356  std::array<double, 4> criterion =
358  if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
359  ++u;
360  ++v;
361  }
362  }
363  if (u == 15 && v == 15) {
364  std::array<double, 4> criterion =
366  if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
367  --u;
368  }
369  }
371  std::array<double, 4> criterion =
373  if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
374  std::pair<double, double> xy1 = hgcalcell_->cellUV2XY1(u, v, placement, 1);
375  std::pair<double, double> xy2 = hgcalcell_->cellUV2XY1(u - 2, v - 1, placement, 1);
376  if (((placement >= HGCalCell::cellPlacementExtra) &&
377  ((((xloc - xy1.first) / (xy2.first - xy1.first)) - ((yloc - xy1.second) / (xy2.second - xy1.second))) >
378  0.0)) ||
379  ((placement < HGCalCell::cellPlacementExtra) &&
380  ((((xloc - xy1.first) / (xy2.first - xy1.first)) - ((yloc - xy1.second) / (xy2.second - xy1.second))) <
381  0.0))) {
382  --u;
383  if ((v - u) >= ncell_[1])
384  --v;
385  } else {
386  --u;
387  --v;
388  v = std::max(v, 0);
389  }
390  }
391  }
392  } else if ((partial != HGCalTypes::WaferFull) && (type == 0)) {
393  if (u == 10 && v == 0) {
394  std::array<double, 4> criterion =
396  if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
397  --u;
398  }
399  }
400  if (u == 10 && v == 21) {
401  std::array<double, 4> criterion =
403  if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
404  --u;
405  --v;
406  }
407  }
410  std::array<double, 4> criterion =
412  if ((criterion[0] * yloc) + (criterion[1] * xloc) < -criterion[2]) {
413  std::pair<double, double> xy1 = hgcalcell_->cellUV2XY1(u, v, placement, 0);
414  std::pair<double, double> xy2 = hgcalcell_->cellUV2XY1(u - 2, v - 1, placement, 0);
415  if (((placement >= HGCalCell::cellPlacementExtra) &&
416  ((((xloc - xy1.first) / (xy2.first - xy1.first)) - ((yloc - xy1.second) / (xy2.second - xy1.second))) >
417  0.0)) ||
418  ((placement < HGCalCell::cellPlacementExtra) &&
419  ((((xloc - xy1.first) / (xy2.first - xy1.first)) - ((yloc - xy1.second) / (xy2.second - xy1.second))) <
420  0.0))) {
421  ++u;
422  ++v;
423  } else {
424  ++u;
425  }
426  }
427  }
428  }
429  if (debug)
430  edm::LogVerbatim("HGCalGeom") << "cellUVFromXY5: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
431  << ":" << v;
432  return std::make_pair(u, v);
433 }
static constexpr std::array< int, 3 > edgeWaferHDBottom
Definition: HGCalTypes.h:88
Log< level::Info, true > LogVerbatim
static constexpr std::array< int, 3 > edgeWaferLDTop
Definition: HGCalTypes.h:81
std::pair< int32_t, int32_t > cellUVFromXY1(double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug) const
Definition: HGCalCellUV.cc:46
static constexpr int32_t WaferLDThree
Definition: HGCalTypes.h:50
static constexpr int32_t WaferFull
Definition: HGCalTypes.h:35
static std::array< double, 4 > maskCut(const int &part, const int &place, const double &waferSize, const double &offset, const bool &v17OrLess)
const double waferSize_
Definition: HGCalCellUV.h:55
static constexpr int32_t WaferHDBottom
Definition: HGCalTypes.h:52
std::unique_ptr< HGCalCell > hgcalcell_
Definition: HGCalCellUV.h:58
#define debug
Definition: HDRShower.cc:19
int32_t ncell_[2]
Definition: HGCalCellUV.h:56
static constexpr int32_t WaferLDTop
Definition: HGCalTypes.h:45
static constexpr int32_t WaferHDTop
Definition: HGCalTypes.h:51
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

◆ cellUVFromXY3()

std::pair< int32_t, int32_t > HGCalCellUV::cellUVFromXY3 ( double  xloc,
double  yloc,
int32_t  placement,
int32_t  type,
bool  extend,
bool  debug 
) const

Definition at line 159 of file HGCalCellUV.cc.

References funct::abs(), reco::ceil(), HGCalCell::cellPlacementExtra, cellX_, cellXTotal_, cellY_, cellYTotal_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), cos60_, cuy::cv, debug, HLT_2024v11_cff::distance, mps_fire::i, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::iv, ncell_, makeMuonMisalignmentScenario::rot, sin60_, sqrt3_, findQualityFiles::v, x, and y.

160  {
161  //--- Using Cube coordinates to find the (u, v)
162  //--- Reverse transform to placement=0, if placement index ≠ 6
163  if (type != 0)
164  type = 1;
165  double cellX = (extend) ? cellXTotal_[type] : cellX_[type];
166  double cellY = (extend) ? cellYTotal_[type] : cellY_[type];
167  double xloc1 = (placement >= HGCalCell::cellPlacementExtra) ? xloc : -xloc;
168  int rot = placement % HGCalCell::cellPlacementExtra;
169  static constexpr std::array<double, 6> fcos = {{1.0, cos60_, -cos60_, -1.0, -cos60_, cos60_}};
170  static constexpr std::array<double, 6> fsin = {{0.0, sin60_, sin60_, 0.0, -sin60_, -sin60_}};
171  double xprime = xloc1 * fcos[rot] - yloc * fsin[rot];
172  double yprime = xloc1 * fsin[rot] + yloc * fcos[rot];
173  double x = xprime + cellX;
174  double y = yprime;
175 
176  x = x / cellX;
177  y = y / cellY;
178 
179  double cu = 2 * x / 3;
180  double cv = -x / 3 + y / 2;
181  double cw = -x / 3 - y / 2;
182 
183  int iu = std::round(cu);
184  int iv = std::round(cv);
185  int iw = std::round(cw);
186 
187  if (iu + iv + iw != 0) {
188  double arr[] = {std::abs(cu - iu), std::abs(cv - iv), std::abs(cw - iw)};
189  int i = std::distance(arr, std::max_element(arr, arr + 3));
190 
191  if (i == 1)
192  iv = (std::round(cv) == std::floor(cv)) ? std::ceil(cv) : std::floor(cv);
193  else if (i == 2)
194  iw = (std::round(cw) == std::floor(cw)) ? std::ceil(cw) : std::floor(cw);
195  }
196 
197  //--- Taking care of extending cells
198  int u(ncell_[type] + iv), v(ncell_[type] - 1 - iw);
199  double xcell = (1.5 * (v - u) + 0.5) * cellX_[type];
200  double ycell = (v + u - 2 * ncell_[type] + 1) * cellY;
201  if (v == -1) {
202  if ((yprime - sqrt3_ * xprime) < (ycell - sqrt3_ * xcell)) {
203  ++v;
204  } else {
205  ++u;
206  ++v;
207  }
208  }
209  if (v - u == ncell_[type]) {
210  if (yprime < ycell) {
211  --v;
212  } else {
213  ++u;
214  }
215  }
216  if (u == 2 * ncell_[type]) {
217  if ((yprime + sqrt3_ * xprime) > (ycell + sqrt3_ * xcell)) {
218  --u;
219  } else {
220  --u;
221  --v;
222  }
223  }
224 
225  if (debug)
226  edm::LogVerbatim("HGCalGeom") << "cellUVFromXY3: Input " << xloc << ":" << yloc << ":" << extend << " Output " << u
227  << ":" << v;
228  return std::make_pair(u, v);
229 }
Log< level::Info, true > LogVerbatim
constexpr int32_t ceil(float num)
static constexpr double sqrt3_
Definition: HGCalCellUV.h:51
double cellXTotal_[2]
Definition: HGCalCellUV.h:57
cv
Definition: cuy.py:363
double cellYTotal_[2]
Definition: HGCalCellUV.h:57
double cellY_[2]
Definition: HGCalCellUV.h:57
double cellX_[2]
Definition: HGCalCellUV.h:57
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr double cos60_
Definition: HGCalCellUV.h:53
#define debug
Definition: HDRShower.cc:19
int32_t ncell_[2]
Definition: HGCalCellUV.h:56
static constexpr double sin60_
Definition: HGCalCellUV.h:52
static constexpr int32_t cellPlacementExtra
Definition: HGCalCell.h:24

◆ cellUVFromXY4() [1/2]

std::pair< int, int > HGCalCellUV::cellUVFromXY4 ( double  xloc,
double  yloc,
int32_t  placement,
int32_t  type,
bool  extend,
bool  debug 
)

Definition at line 231 of file HGCalCellUV.cc.

References cellPosCoarse_, cellPosFine_, cellX_, cellXTotal_, cellY_, debug, and ncell_.

232  {
233  if (type != 0)
234  return cellUVFromXY4(
235  xloc, yloc, ncell_[1], cellX_[1], cellY_[1], cellXTotal_[1], cellY_[1], cellPosCoarse_[placement], extend, debug);
236  else
237  return cellUVFromXY4(
238  xloc, yloc, ncell_[0], cellX_[0], cellY_[0], cellXTotal_[0], cellY_[0], cellPosFine_[placement], extend, debug);
239 }
std::map< std::pair< int32_t, int32_t >, std::pair< double, double > > cellPosFine_[HGCalCell::cellPlacementTotal]
Definition: HGCalCellUV.h:59
double cellXTotal_[2]
Definition: HGCalCellUV.h:57
double cellY_[2]
Definition: HGCalCellUV.h:57
double cellX_[2]
Definition: HGCalCellUV.h:57
#define debug
Definition: HDRShower.cc:19
int32_t ncell_[2]
Definition: HGCalCellUV.h:56
std::map< std::pair< int32_t, int32_t >, std::pair< double, double > > cellPosCoarse_[HGCalCell::cellPlacementTotal]
Definition: HGCalCellUV.h:59
std::pair< int32_t, int32_t > cellUVFromXY4(double xloc, double yloc, int32_t placement, int32_t type, bool extend, bool debug)
Definition: HGCalCellUV.cc:231

◆ cellUVFromXY4() [2/2]

std::pair< int, int > HGCalCellUV::cellUVFromXY4 ( double  xloc,
double  yloc,
int  ncell,
double  cellX,
double  cellY,
double  cellXTotal,
double  cellYTotal,
std::map< std::pair< int, int >, std::pair< double, double > > &  cellPos,
bool  extend,
bool  debug 
)
private

Definition at line 241 of file HGCalCellUV.cc.

References funct::abs(), debug, and sqrt3_.

250  {
251  std::pair<int, int> uv = std::make_pair(-1, -1);
252  std::map<std::pair<int, int>, std::pair<double, double> >::const_iterator itr;
253  for (itr = cellPos.begin(); itr != cellPos.end(); ++itr) {
254  double delX = std::abs(xloc - (itr->second).first);
255  double delY = std::abs(yloc - (itr->second).second);
256  if ((delX < cellX) && (delY < cellY)) {
257  if ((delX < (0.5 * cellX)) || (delY < (2.0 * cellY - sqrt3_ * delX))) {
258  uv = itr->first;
259  break;
260  }
261  }
262  }
263  if ((uv.first < 0) && extend) {
264  for (itr = cellPos.begin(); itr != cellPos.end(); ++itr) {
265  double delX = std::abs(xloc - (itr->second).first);
266  double delY = std::abs(yloc - (itr->second).second);
267  if ((delX < cellXTotal) && (delY < cellYTotal)) {
268  if ((delX < (0.5 * cellXTotal)) || (delY < (2.0 * cellYTotal - sqrt3_ * delX))) {
269  uv = itr->first;
270  break;
271  }
272  }
273  }
274  }
275  if (debug)
276  edm::LogVerbatim("HGCalGeom") << "cellUVFromXY4: Input " << xloc << ":" << yloc << ":" << extend << " Output "
277  << uv.first << ":" << uv.second;
278  return uv;
279 }
Log< level::Info, true > LogVerbatim
static constexpr double sqrt3_
Definition: HGCalCellUV.h:51
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define debug
Definition: HDRShower.cc:19

Member Data Documentation

◆ cellPosCoarse_

std::map<std::pair<int32_t, int32_t>, std::pair<double, double> > HGCalCellUV::cellPosCoarse_[HGCalCell::cellPlacementTotal]
private

Definition at line 59 of file HGCalCellUV.h.

Referenced by cellUVFromXY4(), and HGCalCellUV().

◆ cellPosFine_

std::map<std::pair<int32_t, int32_t>, std::pair<double, double> > HGCalCellUV::cellPosFine_[HGCalCell::cellPlacementTotal]
private

Definition at line 59 of file HGCalCellUV.h.

Referenced by cellUVFromXY4(), and HGCalCellUV().

◆ cellX_

double HGCalCellUV::cellX_[2]
private

Definition at line 57 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY3(), cellUVFromXY4(), and HGCalCellUV().

◆ cellXTotal_

double HGCalCellUV::cellXTotal_[2]
private

Definition at line 57 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY3(), cellUVFromXY4(), and HGCalCellUV().

◆ cellY_

double HGCalCellUV::cellY_[2]
private

Definition at line 57 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY2(), cellUVFromXY3(), cellUVFromXY4(), and HGCalCellUV().

◆ cellYTotal_

double HGCalCellUV::cellYTotal_[2]
private

Definition at line 57 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY2(), cellUVFromXY3(), and HGCalCellUV().

◆ cos60_

constexpr double HGCalCellUV::cos60_ = 0.5
staticprivate

Definition at line 53 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY2(), and cellUVFromXY3().

◆ hgcalcell_

std::unique_ptr<HGCalCell> HGCalCellUV::hgcalcell_
private

Definition at line 58 of file HGCalCellUV.h.

Referenced by cellUVFromXY2(), and HGCalCellUV().

◆ ncell_

int32_t HGCalCellUV::ncell_[2]
private

Definition at line 56 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY2(), cellUVFromXY3(), cellUVFromXY4(), and HGCalCellUV().

◆ sin60_

constexpr double HGCalCellUV::sin60_ = 0.5 * sqrt3_
staticprivate

Definition at line 52 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY2(), and cellUVFromXY3().

◆ sqrt3_

constexpr double HGCalCellUV::sqrt3_ = 1.732050807568877
staticprivate

Definition at line 51 of file HGCalCellUV.h.

Referenced by cellUVFromXY1(), cellUVFromXY2(), cellUVFromXY3(), cellUVFromXY4(), and HGCalCellUV().

◆ waferSize_

const double HGCalCellUV::waferSize_
private

Definition at line 55 of file HGCalCellUV.h.

Referenced by cellUVFromXY2().