1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
|
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
#include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
#include "Geometry/HGCalCommonData/interface/HGCGuardRing.h"
#include <iostream>
//#define EDM_ML_DEBUG
HGCGuardRing::HGCGuardRing(const HGCalDDDConstants& hgc)
: hgcons_(hgc),
modeUV_(hgcons_.geomMode()),
v17OrLess_(hgcons_.v17OrLess()),
waferSize_(hgcons_.waferSize(false)),
sensorSizeOffset_(hgcons_.sensorSizeOffset(false)),
guardRingOffset_(hgcons_.guardRingOffset(false)) {
offset_ = sensorSizeOffset_ + guardRingOffset_;
offsetPartial_ = guardRingOffset_;
xmax_ = 0.5 * waferSize_ - offset_;
ymax_ = xmax_ / sqrt3_;
c22_ = (v17OrLess_) ? HGCalTypes::c22O : HGCalTypes::c22;
c27_ = (v17OrLess_) ? HGCalTypes::c27O : HGCalTypes::c27;
#ifdef EDM_ML_DEBUG
edm::LogVerbatim("HGCSim") << "Creating HGCGuardRing with wafer size " << waferSize_ << ", Offsets "
<< sensorSizeOffset_ << ":" << guardRingOffset_ << ":" << offset_ << ":" << offsetPartial_
<< ", mode " << modeUV_ << ", xmax|ymax " << xmax_ << ":" << ymax_ << " and c22:c27 "
<< c22_ << ":" << c27_;
#endif
}
bool HGCGuardRing::exclude(G4ThreeVector& point, int zside, int frontBack, int layer, int waferU, int waferV) {
bool check(false);
if (hgcons_.waferHexagon8Module()) {
int index = HGCalWaferIndex::waferIndex(layer, waferU, waferV);
int partial = HGCalWaferType::getPartial(index, hgcons_.getParameter()->waferInfoMap_);
#ifdef EDM_ML_DEBUG
edm::LogVerbatim("HGCSim") << "HGCGuardRing::exclude: Layer " << layer << " wafer " << waferU << ":" << waferV
<< " index " << index << " partial " << partial;
#endif
if (partial == HGCalTypes::WaferFull) {
double dx = std::abs(point.x());
double dy = std::abs(point.y());
if (dx > xmax_) {
check = true;
} else if (dy > (2 * ymax_)) {
check = true;
} else {
check = (dx > (sqrt3_ * (2 * ymax_ - dy)));
}
#ifdef EDM_ML_DEBUG
edm::LogVerbatim("HGCSim") << "HGCGuardRing::exclude: Point " << point << " zside " << zside << " layer " << layer
<< " wafer " << waferU << ":" << waferV << " partial type " << partial << ":"
<< HGCalTypes::WaferFull << " x " << dx << ":" << xmax_ << " y " << dy << ":" << ymax_
<< " check " << check;
#endif
} else if (partial > 0) {
int orient = HGCalWaferType::getOrient(index, hgcons_.getParameter()->waferInfoMap_);
#ifdef EDM_ML_DEBUG
edm::LogVerbatim("HGCSim") << "HGCGuardRing::exclude: Orient " << orient << " Mode " << modeUV_;
#endif
if (hgcons_.v16OrLess()) {
std::vector<std::pair<double, double> > wxy =
HGCalWaferMask::waferXY(partial, orient, zside, waferSize_, offset_, 0.0, 0.0, v17OrLess_);
check = !(insidePolygon(point.x(), point.y(), wxy));
#ifdef EDM_ML_DEBUG
std::ostringstream st1;
st1 << "HGCGuardRing::exclude: Point " << point << " Partial/orient/zside/size/offset " << partial << ":"
<< orient << ":" << zside << ":" << waferSize_ << offset_ << " with " << wxy.size() << " points:";
for (unsigned int k = 0; k < wxy.size(); ++k)
st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
edm::LogVerbatim("HGCSim") << st1.str();
#endif
} else {
int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
std::vector<std::pair<double, double> > wxy =
HGCalWaferMask::waferXY(partial, placement, waferSize_, offset_, 0.0, 0.0, v17OrLess_);
check = !(insidePolygon(point.x(), point.y(), wxy));
#ifdef EDM_ML_DEBUG
std::ostringstream st1;
st1 << "HGCGuardRing::exclude: Point " << point << " Partial/frontback/orient/zside/placeemnt/size/offset "
<< partial << ":" << frontBack << ":" << orient << ":" << zside << ":" << placement << ":" << waferSize_
<< offset_ << " with " << wxy.size() << " points:";
for (unsigned int k = 0; k < wxy.size(); ++k)
st1 << " (" << wxy[k].first << ", " << wxy[k].second << ")";
edm::LogVerbatim("HGCSim") << st1.str();
#endif
}
} else {
check = true;
}
}
return check;
}
bool HGCGuardRing::excludePartial(G4ThreeVector& point, int zside, int frontBack, int layer, int waferU, int waferV) {
bool check(false);
if (hgcons_.waferHexagon8Cassette()) {
int index = HGCalWaferIndex::waferIndex(layer, waferU, waferV);
int partial = HGCalWaferType::getPartial(index, hgcons_.getParameter()->waferInfoMap_);
int type = HGCalWaferType::getType(index, hgcons_.getParameter()->waferInfoMap_);
#ifdef EDM_ML_DEBUG
edm::LogVerbatim("HGCSim") << "HGCGuardRing::excludePartial: Layer " << layer << " wafer " << waferU << ":"
<< waferV << " index " << index << " partial " << partial << " type " << type;
#endif
if (partial == HGCalTypes::WaferFull) {
return (check);
} else if (partial < 0) {
return true;
} else {
int orient = HGCalWaferType::getOrient(index, hgcons_.getParameter()->waferInfoMap_);
int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
double dx = point.x();
double dy = point.y();
#ifdef EDM_ML_DEBUG
edm::LogVerbatim("HGCSim") << "HGCGuardRing::excludePartial: zside " << zside << " frontBack " << frontBack
<< " orient " << orient << " placement " << placement << " dx " << dx << " dy " << dy;
#endif
if (type > 0) {
for (int ii = HGCalTypes::WaferPartLDOffset;
ii < (HGCalTypes::WaferPartLDOffset + HGCalTypes::WaferPartLDCount);
ii++) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(ii, placement, waferSize_, offsetPartial_, v17OrLess_);
check |= std::abs(criterion[0] * dy + criterion[1] * dx + criterion[2]) < criterion[3];
}
} else {
for (int ii = HGCalTypes::WaferPartHDOffset;
ii < (HGCalTypes::WaferPartHDOffset + HGCalTypes::WaferPartHDCount);
ii++) {
std::array<double, 4> criterion =
HGCalWaferMask::maskCut(ii, placement, waferSize_, offsetPartial_, v17OrLess_);
check |= std::abs(criterion[0] * dy + criterion[1] * dx + criterion[2]) < criterion[3];
}
}
}
#ifdef EDM_ML_DEBUG
edm::LogVerbatim("HGCSim") << "HGCGuardRing::excludePartial: Point " << point << " zside " << zside << " layer "
<< layer << " wafer " << waferU << ":" << waferV << " partial type " << partial
<< " type " << type << " check " << check;
#endif
}
return check;
}
bool HGCGuardRing::insidePolygon(double x, double y, const std::vector<std::pair<double, double> >& xyv) {
int counter(0);
double x1(xyv[0].first), y1(xyv[0].second);
for (unsigned i1 = 1; i1 <= xyv.size(); i1++) {
unsigned i2 = (i1 % xyv.size());
double x2(xyv[i2].first), y2(xyv[i2].second);
if (y > std::min(y1, y2)) {
if (y <= std::max(y1, y2)) {
if (x <= std::max(x1, x2)) {
if (y1 != y2) {
double xinter = (y - y1) * (x2 - x1) / (y2 - y1) + x1;
if ((x1 == x2) || (x <= xinter))
++counter;
}
}
}
}
x1 = x2;
y1 = y2;
}
if (counter % 2 == 0)
return false;
else
return true;
}
|