CMS 3D CMS Logo

HGCGuardRingPartial.cc
Go to the documentation of this file.
6 #include <iostream>
7 
8 //#define EDM_ML_DEBUG
9 
11  : hgcons_(hgc),
12  modeUV_(hgcons_.geomMode()),
13  v17OrLess_(hgcons_.v17OrLess()),
14  waferSize_(hgcons_.waferSize(false)),
15  guardRingOffset_(hgcons_.getParameter()->guardRingOffset_) {
19 #ifdef EDM_ML_DEBUG
20  edm::LogVerbatim("HGCSim") << "Creating HGCGuardRingPartial with wafer size " << waferSize_ << ", Offsets "
21  << ":" << guardRingOffset_ << ":" << offset_ << ", and mode " << modeUV_
22  << " coefficients " << c22_ << ":" << c27_;
23 #endif
24 }
25 
26 bool HGCGuardRingPartial::exclude(G4ThreeVector& point, int zside, int frontBack, int layer, int waferU, int waferV) {
27  bool check(false);
32  if (partial == HGCalTypes::WaferFull) {
33  return (check);
34  } else {
36  int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
37  double delX = 0.5 * waferSize_;
38  double delY = 2 * delX / sqrt3_;
39  double dx = (zside > 0) ? -point.x() : point.x();
40  double dy = point.y();
41  double tresh = std::abs(offset_ / cos_1[placement]);
42  if (type > 0) {
43  check |= std::abs(dy - (dx * tan_1[placement])) < tresh;
44  check |= std::abs(dy - (dx * tan_1[placement]) + ((HGCalTypes::c10 * delY * 0.5) / cos_1[placement])) < tresh;
45  check |= std::abs(dy * cot_1[placement] - (dx)) < tresh;
46  } else {
47  check |= std::abs((dy * cot_1[placement]) - dx + ((c22_ * delX) / cos_1[placement])) < tresh;
48  check |= std::abs(dy - (dx * tan_1[placement]) - ((c27_ * delY) / cos_1[placement])) < tresh;
49  check |= std::abs(dy - (dx * tan_1[placement]) + ((c27_ * delY) / cos_1[placement])) < tresh;
50  }
51  }
52 #ifdef EDM_ML_DEBUG
53  edm::LogVerbatim("HGCSim") << "HGCGuardRingPartial:: Point " << point << " zside " << zside << " layer " << layer
54  << " wafer " << waferU << ":" << waferV << " partial type " << partial << " type "
55  << type << " check " << check;
56 #endif
57  }
58  return check;
59 }
Log< level::Info, true > LogVerbatim
static constexpr double c27O
Definition: HGCalTypes.h:97
static int32_t cellPlacementIndex(int32_t iz, int32_t frontBack, int32_t orient)
Definition: HGCalCell.cc:237
static int getType(int index, const HGCalParameters::waferInfo_map &wafers)
const double guardRingOffset_
static constexpr double sqrt3_
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 constexpr double c27
Definition: HGCalTypes.h:98
static constexpr double c10
Definition: HGCalTypes.h:107
static int getOrient(int index, const HGCalParameters::waferInfo_map &wafers)
const HGCalGeometryMode::GeometryMode modeUV_
static constexpr int32_t WaferFull
Definition: HGCalTypes.h:35
HGCGuardRingPartial(const HGCalDDDConstants &hgc)
static constexpr double c22O
Definition: HGCalTypes.h:94
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr std::array< double, 12 > cot_1
int32_t waferIndex(int32_t layer, int32_t waferU, int32_t waferV, bool old=false)
static constexpr std::array< double, 12 > tan_1
int32_t waferV(const int32_t index)
static constexpr double c22
Definition: HGCalTypes.h:95
waferInfo_map waferInfoMap_
bool exclude(G4ThreeVector &point, int zside, int frontBack, int layer, int waferU, int waferV)
const HGCalDDDConstants & hgcons_
static constexpr std::array< double, 12 > cos_1
*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