CMS 3D CMS Logo

HGCalNumberingScheme.cc
Go to the documentation of this file.
1 // File: HGCalNumberingScheme.cc
3 // Description: Numbering scheme for High Granularity Calorimeter
15 #include <array>
16 #include <fstream>
17 #include <iostream>
18 
19 //#define EDM_ML_DEBUG
20 
22  const DetId::Detector& det,
23  const std::string& name,
24  const std::string& fileName)
25  : hgcons_(hgc), mode_(hgc.geomMode()), det_(det), name_(name) {
26 #ifdef EDM_ML_DEBUG
27  edm::LogVerbatim("HGCSim") << "Creating HGCalNumberingScheme for " << name_ << " Det " << det_ << " Mode " << mode_
33 #endif
35  if (!fileName.empty()) {
36  edm::FileInPath filetmp1("SimG4CMS/Calo/data/" + fileName);
37  std::string filetmp2 = filetmp1.fullPath();
38  std::ifstream fInput(filetmp2.c_str());
39  if (!fInput.good()) {
40  edm::LogVerbatim("HGCalSim") << "Cannot open file " << filetmp2;
41  } else {
42  char buffer[80];
43  while (fInput.getline(buffer, 80)) {
44  std::vector<std::string> items = CaloSimUtils::splitString(std::string(buffer));
45  if (items.size() == 3) {
46  if (hgcons_.waferHexagon8File()) {
47  int layer = std::atoi(items[0].c_str());
48  int waferU = std::atoi(items[1].c_str());
49  int waferV = std::atoi(items[2].c_str());
50  indices_.emplace_back(HGCalWaferIndex::waferIndex(layer, waferU, waferV, false));
51  } else if (hgcons_.tileTrapezoid()) {
52  int layer = std::atoi(items[0].c_str());
53  int ring = std::atoi(items[1].c_str());
54  int iphi = std::atoi(items[2].c_str());
56  }
57  } else if (items.size() == 1) {
58  int dumpdet = std::atoi(items[0].c_str());
59  dumpDets_.emplace_back(dumpdet);
60  } else if (items.size() == 4) {
61  int idet = std::atoi(items[0].c_str());
62  int layer = std::atoi(items[1].c_str());
63  int zside = std::atoi(items[2].c_str());
64  int cassette = std::atoi(items[3].c_str());
65  dumpCassette_.emplace_back(HGCalCassette::cassetteIndex(idet, layer, zside, cassette));
66  }
67  }
68 #ifdef EDM_ML_DEBUG
69  edm::LogVerbatim("HGCalSim") << "Reads in " << indices_.size() << ":" << dumpDets_.size() << ":"
70  << dumpCassette_.size() << " component information from " << filetmp2
71  << " Layer Offset " << firstLayer_;
72 #endif
73  fInput.close();
74  }
75  }
76 }
77 
79 #ifdef EDM_ML_DEBUG
80  edm::LogVerbatim("HGCSim") << "Deleting HGCalNumberingScheme";
81 #endif
82 }
83 
84 uint32_t HGCalNumberingScheme::getUnitID(int layer, int module, int cell, int iz, const G4ThreeVector& pos, double& wt) {
85  // module is the copy number of the wafer as placed in the layer
86  uint32_t index(0);
87  wt = 1.0;
88 #ifdef EDM_ML_DEBUG
89  edm::LogVerbatim("HGCSim") << "HGCalNumberingScheme:: input Layer " << layer << " Module " << module << " Cell "
90  << cell << " iz " << iz << " Position " << pos;
91 #endif
92  if (hgcons_.waferHexagon8()) {
93  int cellU(0), cellV(0), waferType(-1), waferU(0), waferV(0);
94  bool debug(false);
95  if (cell >= 0) {
99  cellU = HGCalTypes::getUnpackedCellU(cell);
100  cellV = HGCalTypes::getUnpackedCellV(cell);
101  } else if (mode_ != HGCalGeometryMode::Hexagon8) {
102  int zside = (pos.z() > 0) ? 1 : -1;
103  double xx = zside * pos.x();
106  if (!indices_.empty()) {
107  int indx = HGCalWaferIndex::waferIndex(firstLayer_ + layer, wU, wV, false);
108  if (std::find(indices_.begin(), indices_.end(), indx) != indices_.end())
109  debug = true;
110  }
111  if (!dumpDets_.empty()) {
112  if ((std::find(dumpDets_.begin(), dumpDets_.end(), det_) != dumpDets_.end()) &&
114  debug = true;
115  }
116  if (!dumpCassette_.empty()) {
117  int indx = HGCalWaferIndex::waferIndex(firstLayer_ + layer, wU, wV, false);
118  auto ktr = hgcons_.getParameter()->waferInfoMap_.find(indx);
119  if (ktr != hgcons_.getParameter()->waferInfoMap_.end()) {
120  if (std::find(dumpCassette_.begin(),
121  dumpCassette_.end(),
122  HGCalCassette::cassetteIndex(det_, firstLayer_ + layer, zside, (ktr->second).cassette)) !=
123  dumpCassette_.end())
124  debug = true;
125  }
126  }
127  hgcons_.waferFromPosition(xx, pos.y(), zside, layer, waferU, waferV, cellU, cellV, waferType, wt, false, debug);
128  }
129  if (waferType >= 0) {
130  if (hgcons_.waferHexagon8File()) {
131  int type = hgcons_.waferType(layer, waferU, waferV, true);
132  if (type != waferType) {
133 #ifdef EDM_ML_DEBUG
134  edm::LogVerbatim("HGCSim") << "HGCalNumberingScheme:: " << name_ << " Layer|u|v|Index|module|cell " << layer
135  << ":" << waferU << ":" << waferV << ":"
136  << HGCalWaferIndex::waferIndex(layer, waferU, waferV, false) << ":" << module
137  << ":" << cell << " has a type mismatch " << waferType << ":" << type;
138 #endif
140  waferType = type;
141  }
142  }
143  index = HGCSiliconDetId(det_, iz, waferType, layer, waferU, waferV, cellU, cellV).rawId();
144  if (debug) {
145  int zside = (pos.z() > 0) ? 1 : -1;
146  double xx = zside * pos.x();
147  edm::LogVerbatim("HGCSim") << "OK Input " << det_ << ":" << zside << ":" << layer << ":" << xx << ":" << pos.y()
148  << ":" << module << ":" << cell << " WaferType " << waferType << " Wafer " << waferU
149  << ":" << waferV << " Cell " << cellU << ":" << cellV << " wt " << wt << " Mode "
150  << mode_;
151  }
152 #ifdef EDM_ML_DEBUG
153  } else {
154  edm::LogVerbatim("HGCSim") << "Bad WaferType " << waferType << " for Layer:u:v " << layer << ":" << waferU << ":"
155  << waferV;
156 #endif
157  }
158  } else if (hgcons_.tileTrapezoid()) {
159  std::array<int, 3> id = hgcons_.assignCellTrap(pos.x(), pos.y(), pos.z(), layer, false);
160  if (id[2] >= 0) {
161  std::pair<int, int> typm = hgcons_.tileType(layer, id[0], 0);
162  HGCScintillatorDetId detId(id[2], layer, iz * id[0], id[1], false, 0);
163  if (typm.first >= 0) {
164  detId.setType(typm.first);
165  detId.setSiPM(typm.second);
166  }
167  index = detId.rawId();
168  bool debug(false);
169  if (!indices_.empty()) {
170  int indx = HGCalTileIndex::tileIndex(layer, id[0], id[1]);
171  if (std::find(indices_.begin(), indices_.end(), indx) != indices_.end())
172  debug = true;
173  }
174  if (debug)
175  edm::LogVerbatim("HGCSim") << "Radius/Phi " << id[0] << ":" << id[1] << " Type " << id[2] << ":" << typm.first
176  << " SiPM " << typm.second << ":" << hgcons_.tileSiPM(typm.second) << " Layer "
177  << layer << " z " << iz << " " << detId << " wt " << wt << " position " << pos
178  << " R " << pos.perp();
179 #ifdef EDM_ML_DEBUG
180  } else {
181  edm::LogVerbatim("HGCSim") << "Radius/Phi " << id[0] << ":" << id[1] << " Type " << id[2] << " Layer|iz " << layer
182  << ":" << iz << " for i/p Layer " << layer << " module " << module << " cell " << cell
183  << " iz " << iz << " pos " << pos << " wt " << wt << " ERROR";
184 #endif
185  }
186  }
187 #ifdef EDM_ML_DEBUG
188  bool matchOnly = hgcons_.waferHexagon8Module();
190  if (debug)
191  edm::LogVerbatim("HGCSim") << "HGCalNumberingScheme::i/p " << det_ << ":" << layer << ":" << module << ":" << cell
192  << ":" << iz << ":" << pos.x() << ":" << pos.y() << ":" << pos.z() << " ID " << std::hex
193  << index << std::dec << " wt " << wt;
194  bool ok = checkPosition(index, pos, matchOnly, debug);
195  if (matchOnly && (!ok))
196  edm::LogVerbatim("HGCSim") << "HGCalNumberingScheme::i/p " << det_ << ":" << layer << ":" << module << ":" << cell
197  << ":" << iz << ":" << pos.x() << ":" << pos.y() << ":" << pos.z() << " ID " << std::hex
198  << index << std::dec << " wt " << wt << " flag " << ok << " ERROR";
199 #endif
200  return index;
201 }
202 
203 bool HGCalNumberingScheme::checkPosition(uint32_t index, const G4ThreeVector& pos, bool matchOnly, bool debug) const {
204  std::pair<float, float> xy;
205  bool ok(false), iok(true);
206  double z1(0), tolR(14.0), tolZ(1.0);
207  int lay(-1);
208  if (index == 0) {
209  } else if (DetId(index).det() == DetId::HGCalHSi) {
211  lay = id.layer();
213  id.zside(), lay, id.waferU(), id.waferV(), id.cellU(), id.cellV(), false, true, false, false, false);
214  z1 = hgcons_.waferZ(lay, false);
215  ok = true;
216  tolR = 14.0;
217  tolZ = 1.0;
218  } else if (DetId(index).det() == DetId::HGCalHSc) {
220  lay = id.layer();
221  xy = hgcons_.locateCellTrap(id.zside(), lay, id.ietaAbs(), id.iphi(), false, false);
222  z1 = hgcons_.waferZ(lay, false);
223  ok = true;
224  tolR = 50.0;
225  tolZ = 5.0;
226  }
227  if (ok) {
228  double r1 = std::sqrt(xy.first * xy.first + xy.second * xy.second);
229  double r2 = pos.perp();
230  double z2 = std::abs(pos.z());
231  std::pair<double, double> zrange = hgcons_.rangeZ(false);
232  std::pair<double, double> rrange = hgcons_.rangeR(z2, false);
233  bool match = (std::abs(r1 - r2) < tolR) && (std::abs(z1 - z2) < tolZ);
234  bool inok = ((r2 >= rrange.first) && (r2 <= rrange.second) && (z2 >= zrange.first) && (z2 <= zrange.second));
235  bool outok = ((r1 >= rrange.first) && (r1 <= rrange.second) && (z1 >= zrange.first) && (z1 <= zrange.second));
236  std::string ck = (((r1 < rrange.first - tolR) || (r1 > rrange.second + tolR) || (z1 < zrange.first - tolZ) ||
237  (z1 > zrange.second + tolZ))
238  ? "***** ERROR *****"
239  : "");
240  if (matchOnly && match)
241  ck = "";
242  if (!ck.empty())
243  iok = false;
244  if (!(match && inok && outok) || debug) {
245  edm::LogVerbatim("HGCSim") << "HGCalNumberingScheme::Detector " << det_ << " Layer " << lay << " R " << r2 << ":"
246  << r1 << ":" << rrange.first << ":" << rrange.second << " Z " << z2 << ":" << z1 << ":"
247  << zrange.first << ":" << zrange.second << " Match " << match << ":" << inok << ":"
248  << outok << " " << ck;
249  edm::LogVerbatim("HGCSim") << "Original " << pos.x() << ":" << pos.y() << " return " << xy.first << ":"
250  << xy.second;
251  if ((DetId(index).det() == DetId::HGCalEE) || (DetId(index).det() == DetId::HGCalHSi)) {
252  int zside = (pos.z() > 0) ? 1 : -1;
253  double wt(0), xx(zside * pos.x());
254  int waferU, waferV, cellU, cellV, waferType;
255  hgcons_.waferFromPosition(xx, pos.y(), zside, lay, waferU, waferV, cellU, cellV, waferType, wt, false, true);
256  xy = hgcons_.locateCell(zside, lay, waferU, waferV, cellU, cellV, false, true, false, false, true);
257  double dx = (xx - xy.first);
258  double dy = (pos.y() - xy.second);
259  double dR = std::sqrt(dx * dx + dy * dy);
260  if (dR > tolR) {
261  ck = " ***** ERROR *****";
262  iok = false;
263  } else {
264  ck = "";
265  }
266  edm::LogVerbatim("HGCSim") << "HGCalNumberingScheme " << HGCSiliconDetId(index) << " original position " << xx
267  << ":" << pos.y() << " derived " << xy.first << ":" << xy.second << " Difference "
268  << dR << ck;
269  }
270  }
271  }
272  return iok;
273 }
double waferZ(int layer, bool reco) const
Log< level::Info, true > LogVerbatim
static int cassetteIndex(int det, int layer, int zside, int cassette)
void waferFromPosition(const double x, const double y, int &wafer, int &icell, int &celltyp) const
HGCalNumberingScheme()=delete
std::pair< double, double > rangeZ(bool reco) const
const DetId::Detector det_
const HGCalParameters * getParameter() const
static int32_t getUnpackedCellU(int id)
Definition: HGCalTypes.cc:32
std::string fullPath() const
Definition: FileInPath.cc:161
static int32_t getUnpackedU(int id)
Definition: HGCalTypes.cc:16
const HGCalGeometryMode::GeometryMode mode_
int32_t waferU(const int32_t index)
HGCalParameters::waferInfo waferInfo(int lay, int waferU, int waferV) const
const HGCalDDDConstants & hgcons_
int getLayerOffset() const
static int32_t getUnpackedV(int id)
Definition: HGCalTypes.cc:22
bool waferHexagon8() const
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
std::vector< int > indices_
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ietaAbs(uint32_t id)
int32_t tileIndex(int32_t layer, int32_t ring, int32_t phi)
static constexpr int32_t WaferFull
Definition: HGCalTypes.h:35
std::vector< std::string > splitString(const std::string &)
Definition: CaloSimUtils.cc:3
T sqrt(T t)
Definition: SSEVec.h:23
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > dumpDets_
static int32_t getUnpackedType(int id)
Definition: HGCalTypes.cc:14
int tileSiPM(int sipm) const
bool tileTrapezoid() const
std::pair< int, int > tileType(int layer, int ring, int phi) const
std::pair< float, float > locateCellTrap(int zside, int lay, int ieta, int iphi, bool reco, bool debug) const
int32_t waferIndex(int32_t layer, int32_t waferU, int32_t waferV, bool old=false)
bool waferHexagon8Module() const
std::pair< double, double > rangeR(double z, bool reco) const
Definition: DetId.h:17
uint32_t getUnitID(int layer, int module, int cell, int iz, const G4ThreeVector &pos, double &wt)
assigns the det id to a hit
std::array< int, 3 > assignCellTrap(float x, float y, float z, int lay, bool reco) const
#define debug
Definition: HDRShower.cc:19
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
Detector
Definition: DetId.h:24
int waferType(DetId const &id, bool fromFile) const
std::pair< float, float > locateCell(int cell, int lay, int type, bool reco) const
bool waferHexagon8File() const
static int32_t getUnpackedCellV(int id)
Definition: HGCalTypes.cc:34
std::vector< int > dumpCassette_
int32_t waferV(const int32_t index)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
bool checkPosition(uint32_t index, const G4ThreeVector &pos, bool matchOnly, bool debug) const
waferInfo_map waferInfoMap_
const std::string name_