Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-11-22 02:45:25

0001 #include "SimG4CMS/Calo/interface/HGCGuardRingPartial.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
0004 #include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
0005 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0006 #include <iostream>
0007 #include <array>
0008 
0009 //#define EDM_ML_DEBUG
0010 
0011 HGCGuardRingPartial::HGCGuardRingPartial(const HGCalDDDConstants& hgc)
0012     : hgcons_(hgc),
0013       modeUV_(hgcons_.geomMode()),
0014       v17OrLess_(hgcons_.v17OrLess()),
0015       waferSize_(hgcons_.waferSize(false)),
0016       guardRingOffset_(hgcons_.guardRingOffset(false)) {
0017   offset_ = guardRingOffset_;
0018   c22_ = (v17OrLess_) ? HGCalTypes::c22O : HGCalTypes::c22;
0019   c27_ = (v17OrLess_) ? HGCalTypes::c27O : HGCalTypes::c27;
0020 #ifdef EDM_ML_DEBUG
0021   edm::LogVerbatim("HGCSim") << "Creating HGCGuardRingPartial with wafer size " << waferSize_ << ", Offsets "
0022                              << ":" << guardRingOffset_ << ":" << offset_ << ", and mode " << modeUV_
0023                              << " coefficients " << c22_ << ":" << c27_;
0024 #endif
0025 }
0026 
0027 bool HGCGuardRingPartial::exclude(G4ThreeVector& point, int zside, int frontBack, int layer, int waferU, int waferV) {
0028   bool check(false);
0029   if ((modeUV_ == HGCalGeometryMode::Hexagon8Cassette) || (modeUV_ == HGCalGeometryMode::Hexagon8CalibCell)) {
0030     int index = HGCalWaferIndex::waferIndex(layer, waferU, waferV);
0031     int partial = HGCalWaferType::getPartial(index, hgcons_.getParameter()->waferInfoMap_);
0032     int type = HGCalWaferType::getType(index, hgcons_.getParameter()->waferInfoMap_);
0033 #ifdef EDM_ML_DEBUG
0034     edm::LogVerbatim("HGCSim") << "HGCGuardRingPatial:: Layer " << layer << " wafer " << waferU << ":" << waferV
0035                                << " index " << index << " partial " << partial << " type " << type;
0036 #endif
0037     if (partial == HGCalTypes::WaferFull) {
0038       return (check);
0039     } else if (partial < 0) {
0040       return true;
0041     } else {
0042       int orient = HGCalWaferType::getOrient(index, hgcons_.getParameter()->waferInfoMap_);
0043       int placement = HGCalCell::cellPlacementIndex(zside, frontBack, orient);
0044       double dx = point.x();
0045       double dy = point.y();
0046 #ifdef EDM_ML_DEBUG
0047       edm::LogVerbatim("HGCSim") << "HGCGuardRingPatial:: orient " << orient << " placement " << placement << " dx "
0048                                  << dx << " dy " << dy;
0049 #endif
0050       if (type > 0) {
0051         for (int ii = HGCalTypes::WaferPartLDOffset;
0052              ii < (HGCalTypes::WaferPartLDOffset + HGCalTypes::WaferPartLDCount);
0053              ii++) {
0054           std::array<double, 4> criterion = HGCalWaferMask::maskCut(ii, placement, waferSize_, offset_, v17OrLess_);
0055           check |= std::abs(criterion[0] * dy + criterion[1] * dx + criterion[2]) < criterion[3];
0056         }
0057       } else {
0058         for (int ii = HGCalTypes::WaferPartHDOffset;
0059              ii < (HGCalTypes::WaferPartHDOffset + HGCalTypes::WaferPartHDCount);
0060              ii++) {
0061           std::array<double, 4> criterion = HGCalWaferMask::maskCut(ii, placement, waferSize_, offset_, v17OrLess_);
0062           check |= std::abs(criterion[0] * dy + criterion[1] * dx + criterion[2]) < criterion[3];
0063         }
0064       }
0065     }
0066 #ifdef EDM_ML_DEBUG
0067     edm::LogVerbatim("HGCSim") << "HGCGuardRingPartial:: Point " << point << " zside " << zside << " layer " << layer
0068                                << " wafer " << waferU << ":" << waferV << " partial type " << partial << " type "
0069                                << type << " check " << check;
0070 #endif
0071   }
0072   return check;
0073 }