CMS 3D CMS Logo

HGCGuardRing.cc
Go to the documentation of this file.
5 #include <iostream>
6 
7 //#define EDM_ML_DEBUG
8 
10  : hgcons_(hgc),
11  modeUV_(hgcons_.geomMode()),
12  v17OrLess_(hgcons_.v17OrLess()),
13  waferSize_(hgcons_.waferSize(false)),
14  sensorSizeOffset_(hgcons_.getParameter()->sensorSizeOffset_),
15  guardRingOffset_(hgcons_.getParameter()->guardRingOffset_) {
17  xmax_ = 0.5 * (waferSize_ - offset_);
18  ymax_ = xmax_ / sqrt3_;
19 #ifdef EDM_ML_DEBUG
20  edm::LogVerbatim("HGCSim") << "Creating HGCGuardRing with wafer size " << waferSize_ << ", Offsets "
21  << sensorSizeOffset_ << ":" << guardRingOffset_ << ":" << offset_ << ", and mode "
22  << modeUV_ << " xmax|ymax " << xmax_ << ":" << ymax_;
23 #endif
24 }
25 
26 bool HGCGuardRing::exclude(G4ThreeVector& point, int zside, int frontBack, int layer, int waferU, int waferV) {
27  bool check(false);
31 #ifdef EDM_ML_DEBUG
32  edm::LogVerbatim("HGCSim") << "HGCGuardRing:: Layer " << layer << " wafer " << waferU << ":" << waferV << " index "
33  << index << " partial " << partial;
34 #endif
35  if (partial == HGCalTypes::WaferFull) {
36  double dx = std::abs(point.x());
37  double dy = std::abs(point.y());
38  if (dx > xmax_) {
39  check = true;
40  } else if (dy > (2 * ymax_)) {
41  check = true;
42  } else {
43  check = (dx > (sqrt3_ * (2 * ymax_ - dy)));
44  }
45 #ifdef EDM_ML_DEBUG
46  edm::LogVerbatim("HGCSim") << "HGCGuardRing:: Point " << point << " zside " << zside << " layer " << layer
47  << " wafer " << waferU << ":" << waferV << " partial type " << partial << ":"
48  << HGCalTypes::WaferFull << " x " << dx << ":" << xmax_ << " y " << dy << ":" << ymax_
49  << " check " << check;
50 #endif
51  } else if (partial > 0) {
53 #ifdef EDM_ML_DEBUG
54  edm::LogVerbatim("HGCSim") << "HGCGuardRing:: Orient " << orient << " Mode " << modeUV_;
55 #endif
57  std::vector<std::pair<double, double> > wxy =
59  check = !(insidePolygon(point.x(), point.y(), wxy));
60 #ifdef EDM_ML_DEBUG
61  std::ostringstream st1;
62  st1 << "HGCGuardRing:: Point " << point << " Partial/orient/zside/size/offset " << partial << ":" << orient
63  << ":" << zside << ":" << waferSize_ << offset_ << " with " << wxy.size() << " points:";
64  for (unsigned int k = 0; k < wxy.size(); ++k)
65  st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
66  edm::LogVerbatim("HGCSim") << st1.str();
67 #endif
68  } else {
69  int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
70  std::vector<std::pair<double, double> > wxy =
71  HGCalWaferMask::waferXY(partial, placement, waferSize_, offset_, 0.0, 0.0, v17OrLess_);
72  check = !(insidePolygon(point.x(), point.y(), wxy));
73 #ifdef EDM_ML_DEBUG
74  std::ostringstream st1;
75  st1 << "HGCGuardRing:: Point " << point << " Partial/frontback/orient/zside/placeemnt/size/offset " << partial
76  << ":" << frontBack << ":" << orient << ":" << zside << ":" << placement << ":" << waferSize_ << offset_
77  << " with " << wxy.size() << " points:";
78  for (unsigned int k = 0; k < wxy.size(); ++k)
79  st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
80  edm::LogVerbatim("HGCSim") << st1.str();
81 #endif
82  }
83  } else {
84  check = true;
85  }
86  }
87  return check;
88 }
89 
90 bool HGCGuardRing::insidePolygon(double x, double y, const std::vector<std::pair<double, double> >& xyv) {
91  int counter(0);
92  double x1(xyv[0].first), y1(xyv[0].second);
93  for (unsigned i1 = 1; i1 <= xyv.size(); i1++) {
94  unsigned i2 = (i1 % xyv.size());
95  double x2(xyv[i2].first), y2(xyv[i2].second);
96  if (y > std::min(y1, y2)) {
97  if (y <= std::max(y1, y2)) {
98  if (x <= std::max(x1, x2)) {
99  if (y1 != y2) {
100  double xinter = (y - y1) * (x2 - x1) / (y2 - y1) + x1;
101  if ((x1 == x2) || (x <= xinter))
102  ++counter;
103  }
104  }
105  }
106  }
107  x1 = x2;
108  y1 = y2;
109  }
110 
111  if (counter % 2 == 0)
112  return false;
113  else
114  return true;
115 }
Log< level::Info, true > LogVerbatim
static int32_t cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient)
Definition: HGCalCell.cc:237
const double sensorSizeOffset_
Definition: HGCGuardRing.h:20
double ymax_
Definition: HGCGuardRing.h:21
const HGCalParameters * getParameter() const
int32_t waferU(const int32_t index)
static int getPartial(int index, const HGCalParameters::waferInfo_map &wafers)
int zside(DetId const &)
static std::vector< std::pair< double, double > > waferXY(const int &part, const int &orient, const int &zside, const double &waferSize, const double &offset, const double &xpos, const double &ypos, const bool &v17)
static int getOrient(int index, const HGCalParameters::waferInfo_map &wafers)
U second(std::pair< T, U > const &p)
const double waferSize_
Definition: HGCGuardRing.h:20
static constexpr int32_t WaferFull
Definition: HGCalTypes.h:35
const double guardRingOffset_
Definition: HGCGuardRing.h:20
bool exclude(G4ThreeVector &point, int zside, int frontBack, int layer, int waferU, int waferV)
Definition: HGCGuardRing.cc:26
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double xmax_
Definition: HGCGuardRing.h:21
int32_t waferIndex(int32_t layer, int32_t waferU, int32_t waferV, bool old=false)
bool waferHexagon8Module() const
const HGCalDDDConstants & hgcons_
Definition: HGCGuardRing.h:17
static bool insidePolygon(double x, double y, const std::vector< std::pair< double, double > > &xyv)
Definition: HGCGuardRing.cc:90
HGCGuardRing(const HGCalDDDConstants &hgc)
Definition: HGCGuardRing.cc:9
const bool v17OrLess_
Definition: HGCGuardRing.h:19
double offset_
Definition: HGCGuardRing.h:21
static std::atomic< unsigned int > counter
const HGCalGeometryMode::GeometryMode modeUV_
Definition: HGCGuardRing.h:18
int32_t waferV(const int32_t index)
float x
static constexpr double sqrt3_
Definition: HGCGuardRing.h:16
waferInfo_map waferInfoMap_
*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