CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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  if (partial == HGCalTypes::WaferFull) {
32  double dx = std::abs(point.x());
33  double dy = std::abs(point.y());
34  if (dx > xmax_) {
35  check = true;
36  } else if (dy > (2 * ymax_)) {
37  check = true;
38  } else {
39  check = (dx > (sqrt3_ * (2 * ymax_ - dy)));
40  }
41 #ifdef EDM_ML_DEBUG
42  edm::LogVerbatim("HGCSim") << "HGCGuardRing:: Point " << point << " zside " << zside << " layer " << layer
43  << " wafer " << waferU << ":" << waferV << " partial type " << partial << ":"
44  << HGCalTypes::WaferFull << " x " << dx << ":" << xmax_ << " y " << dy << ":" << ymax_
45  << " check " << check;
46 #endif
47  } else {
50  std::vector<std::pair<double, double> > wxy =
52  check = !(insidePolygon(point.x(), point.y(), wxy));
53 #ifdef EDM_ML_DEBUG
54  std::ostringstream st1;
55  st1 << "HGCGuardRing:: Point " << point << " Partial/orient/zside/size/offset " << partial << ":" << orient
56  << ":" << zside << ":" << waferSize_ << offset_ << " with " << wxy.size() << " points:";
57  for (unsigned int k = 0; k < wxy.size(); ++k)
58  st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
59  edm::LogVerbatim("HGCSim") << st1.str();
60 #endif
61  } else {
62  int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
63  std::vector<std::pair<double, double> > wxy =
64  HGCalWaferMask::waferXY(partial, placement, waferSize_, offset_, 0.0, 0.0, v17OrLess_);
65  check = !(insidePolygon(point.x(), point.y(), wxy));
66 #ifdef EDM_ML_DEBUG
67  std::ostringstream st1;
68  st1 << "HGCGuardRing:: Point " << point << " Partial/frontback/orient/zside/placeemnt/size/offset " << partial
69  << ":" << frontBack << ":" << orient << ":" << zside << ":" << placement << ":" << waferSize_ << offset_
70  << " with " << wxy.size() << " points:";
71  for (unsigned int k = 0; k < wxy.size(); ++k)
72  st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
73  edm::LogVerbatim("HGCSim") << st1.str();
74 #endif
75  }
76  }
77  }
78  return check;
79 }
80 
81 bool HGCGuardRing::insidePolygon(double x, double y, const std::vector<std::pair<double, double> >& xyv) {
82  int counter(0);
83  double x1(xyv[0].first), y1(xyv[0].second);
84  for (unsigned i1 = 1; i1 <= xyv.size(); i1++) {
85  unsigned i2 = (i1 % xyv.size());
86  double x2(xyv[i2].first), y2(xyv[i2].second);
87  if (y > std::min(y1, y2)) {
88  if (y <= std::max(y1, y2)) {
89  if (x <= std::max(x1, x2)) {
90  if (y1 != y2) {
91  double xinter = (y - y1) * (x2 - x1) / (y2 - y1) + x1;
92  if ((x1 == x2) || (x <= xinter))
93  ++counter;
94  }
95  }
96  }
97  }
98  x1 = x2;
99  y1 = y2;
100  }
101 
102  if (counter % 2 == 0)
103  return false;
104  else
105  return true;
106 }
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)
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:81
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