File indexing completed on 2025-06-03 00:12:04
0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
0003 #include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
0004 #include "Geometry/HGCalCommonData/interface/HGCGuardRing.h"
0005 #include <iostream>
0006
0007
0008
0009 HGCGuardRing::HGCGuardRing(const HGCalDDDConstants& hgc)
0010 : hgcons_(hgc),
0011 modeUV_(hgcons_.geomMode()),
0012 v17OrLess_(hgcons_.v17OrLess()),
0013 waferSize_(hgcons_.waferSize(false)),
0014 sensorSizeOffset_(hgcons_.sensorSizeOffset(false)),
0015 guardRingOffset_(hgcons_.guardRingOffset(false)) {
0016 offset_ = sensorSizeOffset_ + guardRingOffset_;
0017 offsetPartial_ = guardRingOffset_;
0018 xmax_ = 0.5 * waferSize_ - offset_;
0019 ymax_ = xmax_ / sqrt3_;
0020 c22_ = (v17OrLess_) ? HGCalTypes::c22O : HGCalTypes::c22;
0021 c27_ = (v17OrLess_) ? HGCalTypes::c27O : HGCalTypes::c27;
0022 #ifdef EDM_ML_DEBUG
0023 edm::LogVerbatim("HGCSim") << "Creating HGCGuardRing with wafer size " << waferSize_ << ", Offsets "
0024 << sensorSizeOffset_ << ":" << guardRingOffset_ << ":" << offset_ << ":" << offsetPartial_
0025 << ", mode " << modeUV_ << ", xmax|ymax " << xmax_ << ":" << ymax_ << " and c22:c27 "
0026 << c22_ << ":" << c27_;
0027 #endif
0028 }
0029
0030 bool HGCGuardRing::exclude(G4ThreeVector& point, int zside, int frontBack, int layer, int waferU, int waferV) {
0031 bool check(false);
0032 if (hgcons_.waferHexagon8Module()) {
0033 int index = HGCalWaferIndex::waferIndex(layer, waferU, waferV);
0034 int partial = HGCalWaferType::getPartial(index, hgcons_.getParameter()->waferInfoMap_);
0035 #ifdef EDM_ML_DEBUG
0036 edm::LogVerbatim("HGCSim") << "HGCGuardRing::exclude: Layer " << layer << " wafer " << waferU << ":" << waferV
0037 << " index " << index << " partial " << partial;
0038 #endif
0039 if (partial == HGCalTypes::WaferFull) {
0040 double dx = std::abs(point.x());
0041 double dy = std::abs(point.y());
0042 if (dx > xmax_) {
0043 check = true;
0044 } else if (dy > (2 * ymax_)) {
0045 check = true;
0046 } else {
0047 check = (dx > (sqrt3_ * (2 * ymax_ - dy)));
0048 }
0049 #ifdef EDM_ML_DEBUG
0050 edm::LogVerbatim("HGCSim") << "HGCGuardRing::exclude: Point " << point << " zside " << zside << " layer " << layer
0051 << " wafer " << waferU << ":" << waferV << " partial type " << partial << ":"
0052 << HGCalTypes::WaferFull << " x " << dx << ":" << xmax_ << " y " << dy << ":" << ymax_
0053 << " check " << check;
0054 #endif
0055 } else if (partial > 0) {
0056 int orient = HGCalWaferType::getOrient(index, hgcons_.getParameter()->waferInfoMap_);
0057 #ifdef EDM_ML_DEBUG
0058 edm::LogVerbatim("HGCSim") << "HGCGuardRing::exclude: Orient " << orient << " Mode " << modeUV_;
0059 #endif
0060 if (hgcons_.v16OrLess()) {
0061 std::vector<std::pair<double, double> > wxy =
0062 HGCalWaferMask::waferXY(partial, orient, zside, waferSize_, offset_, 0.0, 0.0, v17OrLess_);
0063 check = !(insidePolygon(point.x(), point.y(), wxy));
0064 #ifdef EDM_ML_DEBUG
0065 std::ostringstream st1;
0066 st1 << "HGCGuardRing::exclude: Point " << point << " Partial/orient/zside/size/offset " << partial << ":"
0067 << orient << ":" << zside << ":" << waferSize_ << offset_ << " with " << wxy.size() << " points:";
0068 for (unsigned int k = 0; k < wxy.size(); ++k)
0069 st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
0070 edm::LogVerbatim("HGCSim") << st1.str();
0071 #endif
0072 } else {
0073 int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
0074 std::vector<std::pair<double, double> > wxy =
0075 HGCalWaferMask::waferXY(partial, placement, waferSize_, offset_, 0.0, 0.0, v17OrLess_);
0076 check = !(insidePolygon(point.x(), point.y(), wxy));
0077 #ifdef EDM_ML_DEBUG
0078 std::ostringstream st1;
0079 st1 << "HGCGuardRing::exclude: Point " << point << " Partial/frontback/orient/zside/placeemnt/size/offset "
0080 << partial << ":" << frontBack << ":" << orient << ":" << zside << ":" << placement << ":" << waferSize_
0081 << offset_ << " with " << wxy.size() << " points:";
0082 for (unsigned int k = 0; k < wxy.size(); ++k)
0083 st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
0084 edm::LogVerbatim("HGCSim") << st1.str();
0085 #endif
0086 }
0087 } else {
0088 check = true;
0089 }
0090 }
0091 return check;
0092 }
0093
0094 bool HGCGuardRing::excludePartial(G4ThreeVector& point, int zside, int frontBack, int layer, int waferU, int waferV) {
0095 bool check(false);
0096 if (hgcons_.waferHexagon8Cassette()) {
0097 int index = HGCalWaferIndex::waferIndex(layer, waferU, waferV);
0098 int partial = HGCalWaferType::getPartial(index, hgcons_.getParameter()->waferInfoMap_);
0099 int type = HGCalWaferType::getType(index, hgcons_.getParameter()->waferInfoMap_);
0100 #ifdef EDM_ML_DEBUG
0101 edm::LogVerbatim("HGCSim") << "HGCGuardRing::excludePartial: Layer " << layer << " wafer " << waferU << ":"
0102 << waferV << " index " << index << " partial " << partial << " type " << type;
0103 #endif
0104 if (partial == HGCalTypes::WaferFull) {
0105 return (check);
0106 } else if (partial < 0) {
0107 return true;
0108 } else {
0109 int orient = HGCalWaferType::getOrient(index, hgcons_.getParameter()->waferInfoMap_);
0110 int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
0111 double dx = point.x();
0112 double dy = point.y();
0113 #ifdef EDM_ML_DEBUG
0114 edm::LogVerbatim("HGCSim") << "HGCGuardRing::excludePartial: zside " << zside << " frontBack " << frontBack
0115 << " orient " << orient << " placement " << placement << " dx " << dx << " dy " << dy;
0116 #endif
0117 if (type > 0) {
0118 for (int ii = HGCalTypes::WaferPartLDOffset;
0119 ii < (HGCalTypes::WaferPartLDOffset + HGCalTypes::WaferPartLDCount);
0120 ii++) {
0121 std::array<double, 4> criterion =
0122 HGCalWaferMask::maskCut(ii, placement, waferSize_, offsetPartial_, v17OrLess_);
0123 check |= std::abs(criterion[0] * dy + criterion[1] * dx + criterion[2]) < criterion[3];
0124 }
0125 } else {
0126 for (int ii = HGCalTypes::WaferPartHDOffset;
0127 ii < (HGCalTypes::WaferPartHDOffset + HGCalTypes::WaferPartHDCount);
0128 ii++) {
0129 std::array<double, 4> criterion =
0130 HGCalWaferMask::maskCut(ii, placement, waferSize_, offsetPartial_, v17OrLess_);
0131 check |= std::abs(criterion[0] * dy + criterion[1] * dx + criterion[2]) < criterion[3];
0132 }
0133 }
0134 }
0135 #ifdef EDM_ML_DEBUG
0136 edm::LogVerbatim("HGCSim") << "HGCGuardRing::excludePartial: Point " << point << " zside " << zside << " layer "
0137 << layer << " wafer " << waferU << ":" << waferV << " partial type " << partial
0138 << " type " << type << " check " << check;
0139 #endif
0140 }
0141 return check;
0142 }
0143
0144 bool HGCGuardRing::insidePolygon(double x, double y, const std::vector<std::pair<double, double> >& xyv) {
0145 int counter(0);
0146 double x1(xyv[0].first), y1(xyv[0].second);
0147 for (unsigned i1 = 1; i1 <= xyv.size(); i1++) {
0148 unsigned i2 = (i1 % xyv.size());
0149 double x2(xyv[i2].first), y2(xyv[i2].second);
0150 if (y > std::min(y1, y2)) {
0151 if (y <= std::max(y1, y2)) {
0152 if (x <= std::max(x1, x2)) {
0153 if (y1 != y2) {
0154 double xinter = (y - y1) * (x2 - x1) / (y2 - y1) + x1;
0155 if ((x1 == x2) || (x <= xinter))
0156 ++counter;
0157 }
0158 }
0159 }
0160 }
0161 x1 = x2;
0162 y1 = y2;
0163 }
0164
0165 if (counter % 2 == 0)
0166 return false;
0167 else
0168 return true;
0169 }