CMS 3D CMS Logo

HGCalTestGuardRing.cc
Go to the documentation of this file.
3 
7 
14 
17 
22 
31 
32 #include <fstream>
33 #include <sstream>
34 #include <string>
35 #include <map>
36 
38 public:
40  ~HGCalTestGuardRing() override = default;
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
42 
43 protected:
44  void analyze(edm::Event const&, edm::EventSetup const&) override;
45  void beginJob() override {}
46  void endJob() override {}
47 
48 private:
50  const double guardRingOffset_;
52  std::map<HGCSiliconDetId, int> waferID_;
53 };
54 
56  : nameSense_(ps.getParameter<std::string>("nameSense")),
57  waferFile_(ps.getParameter<std::string>("waferFile")),
58  guardRingOffset_(ps.getParameter<double>("guardRingOffset")),
60  DetId::Detector det = (nameSense_ != "HGCalHESiliconSensitive") ? DetId::HGCalEE : DetId::HGCalHSi;
61  edm::LogVerbatim("HGCalSim") << "Test Guard Ring Offset " << guardRingOffset_ << " for " << nameSense_ << ":" << det
62  << " for wafers read from file " << waferFile_;
63  if (!waferFile_.empty()) {
64  std::string thick[4] = {"h120", "l200", "l300", "h200"};
66  const int partTypeH[6] = {HGCalTypes::WaferFull,
72  const int partTypeL[7] = {HGCalTypes::WaferFull,
79  edm::FileInPath filetmp("SimG4CMS/Calo/data/" + waferFile_);
80  std::string fileName = filetmp.fullPath();
81  std::ifstream fInput(fileName.c_str());
82  if (!fInput.good()) {
83  edm::LogVerbatim("HGCalSim") << "Cannot open file " << fileName;
84  } else {
85  char buffer[80];
86  while (fInput.getline(buffer, 80)) {
87  std::vector<std::string> items = CaloSimUtils::splitString(std::string(buffer));
88  if (items.size() > 6) {
89  int layer = std::atoi(items[0].c_str());
90  int waferU = std::atoi(items[4].c_str());
91  int waferV = std::atoi(items[5].c_str());
92  int thck = static_cast<int>(std::find(thick, thick + 4, items[2]) - thick);
93  int type = (thck < 4) ? addType[thck] : 0;
94  HGCSiliconDetId id(det, -1, type, layer, waferU, waferV, 0, 0);
95  int orient = std::atoi(items[5].c_str());
96  int part = std::atoi(items[1].c_str());
97  if (part >= 0) {
99  part = partTypeH[part];
100  else
101  part = partTypeL[part];
102  }
103  waferID_[id] = orient * 100 + part;
104 #ifdef EDM_ML_DEBUG
105  edm::LogVerbatim("HGCalSim") << "HGCalTestGuardRing::Reads " << id << " Orientation:Partial " << orient << ":"
106  << part;
107 #endif
108  }
109  }
110  edm::LogVerbatim("HGCalSim") << "HGCalTestGuardRing::Reads in " << waferID_.size() << " wafer information from "
111  << fileName;
112  fInput.close();
113  }
114  }
115 }
116 
119  desc.add<std::string>("nameSense", "HGCalEESensitive");
120  desc.add<std::string>("waferFile", "testWafersEE.txt");
121  desc.add<double>("guardRingOffset", 1.0);
122  descriptions.add("hgcalTestGuardRingEE", desc);
123 }
124 
126  // get HGCalGeometry
127  const HGCalGeometry* geom = &iS.getData(geomToken_);
128  const HGCalDDDConstants& hgc = geom->topology().dddConstants();
129  double waferSize = hgc.waferSize(false);
130  HGCalCell wafer(waferSize, hgc.getUVMax(0), hgc.getUVMax(1));
131  const bool v17OrLess = hgc.v17OrLess();
132 
133  // get the hit collection
134  edm::LogVerbatim("HGCalSim") << "HGCalTestGuardRing: Wafer Szie " << waferSize << " v17OrLess " << v17OrLess;
135 
136  // Loop over all IDs
137  int all(0), allSi(0), good(0);
138  for (std::map<HGCSiliconDetId, int>::const_iterator itr = waferID_.begin(); itr != waferID_.end(); ++itr) {
139  HGCSiliconDetId id = itr->first;
140  ++all;
141  if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
142  ++allSi;
143  if (((id.det() == DetId::HGCalEE) && (nameSense_ == "HGCalEESensitive")) ||
144  ((id.det() == DetId::HGCalHSi) && (nameSense_ == "HGCalHESiliconSensitive"))) {
145  int partial = ((itr->second) % 100);
146  int orient = (((itr->second) / 100) % 100);
147  int type = id.type();
148  int nCells = hgc.getUVMax(type);
149  for (int u = 0; u < 2 * nCells; ++u) {
150  for (int v = 0; v < 2 * nCells; ++v) {
151  if (((v - u) < nCells) && ((u - v) <= nCells)) {
152  HGCSiliconDetId hid(id.det(), id.zside(), id.type(), id.layer(), id.waferU(), id.waferV(), u, v);
153  bool valid = (geom->topology()).valid(static_cast<DetId>(hid));
154  if (valid) {
155  ++good;
156  int placeIndex = wafer.cellPlacementIndex(1, HGCalTypes::waferFrontBack(0), orient);
157  std::pair<double, double> xy = wafer.cellUV2XY1(u, v, placeIndex, type);
158  std::vector<std::pair<double, double> > wxy1 =
159  HGCalWaferMask::waferXY(partial, orient, -1, waferSize, 0.0, 0.0, 0.0, v17OrLess);
160  bool check1 = HGCGuardRing::insidePolygon(xy.first, xy.second, wxy1);
161  std::ostringstream st1;
162  for (unsigned int k1 = 0; k1 < wxy1.size(); ++k1)
163  st1 << " (" << wxy1[k1].first << ", " << wxy1[k1].second << ")";
164  edm::LogVerbatim("HGCSim")
165  << "First " << hid << " Type:Partial:Orient:Place " << type << ":" << partial << ":" << orient
166  << ":" << placeIndex << " Boundary with " << wxy1.size() << " points: " << st1.str() << " check "
167  << check1 << " for (" << xy.first << ", " << xy.second << ")";
168 
169  std::vector<std::pair<double, double> > wxy2 =
170  HGCalWaferMask::waferXY(partial, orient, -1, waferSize, guardRingOffset_, 0.0, 0.0, v17OrLess);
171  bool check2 = HGCGuardRing::insidePolygon(xy.first, xy.second, wxy2);
172  std::ostringstream st2;
173  for (unsigned int k1 = 0; k1 < wxy2.size(); ++k1)
174  st2 << " (" << wxy2[k1].first << ", " << wxy2[k1].second << ")";
175  edm::LogVerbatim("HGCSim") << "Second Offset " << guardRingOffset_ << " Boundary with " << wxy2.size()
176  << " points: " << st2.str() << " check " << check2 << " for (" << xy.first
177  << ", " << xy.second << ")";
178  }
179  }
180  }
181  }
182  }
183  }
184  }
185  edm::LogVerbatim("HGCalSim") << "Total hits = " << all << " Good Silicon DetIds = " << allSi << ":" << good;
186 }
187 
188 //define this as a plug-in
Log< level::Info, true > LogVerbatim
static constexpr int32_t WaferHalf2
Definition: HGCalTypes.h:43
static constexpr int32_t WaferFive2
Definition: HGCalTypes.h:44
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
void beginJob() override
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static constexpr int32_t WaferHD200
Definition: HGCalTypes.h:33
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
static constexpr int32_t waferFrontBack(int32_t index)
Definition: HGCalTypes.h:138
int32_t waferU(const int32_t index)
std::map< HGCSiliconDetId, int > waferID_
const double guardRingOffset_
static constexpr int32_t WaferThree
Definition: HGCalTypes.h:42
static constexpr int32_t WaferLD300
Definition: HGCalTypes.h:32
int zside(DetId const &)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
static constexpr int32_t WaferSemi2
Definition: HGCalTypes.h:41
const edm::ESGetToken< HGCalGeometry, IdealGeometryRecord > geomToken_
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)
U second(std::pair< T, U > const &p)
static constexpr int32_t WaferFull
Definition: HGCalTypes.h:35
HGCalTestGuardRing(const edm::ParameterSet &ps)
static constexpr int32_t WaferHalf
Definition: HGCalTypes.h:39
std::vector< std::string > splitString(const std::string &)
Definition: CaloSimUtils.cc:3
void endJob() override
void analyze(edm::Event const &, edm::EventSetup const &) override
static constexpr int32_t WaferHD120
Definition: HGCalTypes.h:30
~HGCalTestGuardRing() override=default
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const std::string nameSense_
const std::string waferFile_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static constexpr int32_t WaferChopTwoM
Definition: HGCalTypes.h:38
part
Definition: HCALResponse.h:20
Detector
Definition: DetId.h:24
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static bool insidePolygon(double x, double y, const std::vector< std::pair< double, double > > &xyv)
Definition: HGCGuardRing.cc:90
HLT enums.
ALPAKA_FN_ACC ALPAKA_FN_INLINE void uint32_t const uint32_t CACellT< TrackerTraits > uint32_t * nCells
int32_t waferV(const int32_t index)
static constexpr int32_t WaferFive
Definition: HGCalTypes.h:36
static constexpr int32_t WaferSemi
Definition: HGCalTypes.h:40
static constexpr int32_t WaferLD200
Definition: HGCalTypes.h:31