CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Attributes
BTLNumberingScheme Class Reference

#include <BTLNumberingScheme.h>

Inheritance diagram for BTLNumberingScheme:
MTDNumberingScheme

Public Member Functions

 BTLNumberingScheme ()
 
uint32_t getUnitID (const MTDBaseNumber &baseNumber) const override
 
 ~BTLNumberingScheme () override
 
- Public Member Functions inherited from MTDNumberingScheme
 MTDNumberingScheme ()
 
virtual ~MTDNumberingScheme ()
 

Static Public Attributes

static constexpr std::array< uint32_t, BTLDetId::kRUPerTypeV2 *BTLDetId::kCrystalTypesglobalru2ru
 
static constexpr std::array< uint32_t, BTLDetId::kRUPerTypeV2 *BTLDetId::kCrystalTypesglobalru2type
 
static constexpr uint32_t kBTLcrystalLevel = 9
 
static constexpr uint32_t kBTLmoduleLevel = 8
 
static constexpr std::array< uint32_t, BTLDetId::kModulesPerRUV2negModCopy
 

Detailed Description

Definition at line 7 of file BTLNumberingScheme.h.

Constructor & Destructor Documentation

◆ BTLNumberingScheme()

BTLNumberingScheme::BTLNumberingScheme ( )

Definition at line 8 of file BTLNumberingScheme.cc.

References LogDebug.

9  LogDebug("MTDGeom") << "Creating BTLNumberingScheme";
10 }
#define LogDebug(id)

◆ ~BTLNumberingScheme()

BTLNumberingScheme::~BTLNumberingScheme ( )
override

Definition at line 12 of file BTLNumberingScheme.cc.

References LogDebug.

12 { LogDebug("MTDGeom") << "Deleting BTLNumberingScheme"; }
#define LogDebug(id)

Member Function Documentation

◆ getUnitID()

uint32_t BTLNumberingScheme::getUnitID ( const MTDBaseNumber baseNumber) const
overridevirtual

Implements MTDNumberingScheme.

Definition at line 14 of file BTLNumberingScheme.cc.

References BTLDetId::geographicalId(), MTDBaseNumber::getCopyNumber(), MTDBaseNumber::getLevelName(), MTDBaseNumber::getLevels(), globalru2ru, globalru2type, BTLDetId::HALF_ROD, kBTLcrystalLevel, kBTLmoduleLevel, BTLDetId::kCrystalsPerModuleV2, BTLDetId::kModulesPerRUV2, BTLDetId::kRUPerTypeV2, LogDebug, mergeVDriftHistosByStation::name, negModCopy, DetId::rawId(), BTLDetId::v2, and ecaldqm::zside().

Referenced by PrintMTDSens::dumpTouch().

14  {
15  uint32_t intindex(0);
16  const uint32_t nLevels(baseNumber.getLevels());
17 
18  LogDebug("MTDGeom") << "BTLNumberingScheme geometry levels = " << nLevels;
19 
20  uint32_t zside(999), rodCopy(0), runitCopy(0), modCopy(0), modtyp(0), crystal(0);
21 
22  bool isDD4hepOK(false);
23  if (nLevels == kBTLcrystalLevel + 1) {
24  if (baseNumber.getLevelName(9) == "world_volume_1") {
25  isDD4hepOK = true;
26  }
27  }
28 
29 #ifdef EDM_ML_DEBUG
30  LogDebug("MTDGeom") << "BTLNumberingScheme::getUnitID(): isDD4hep " << isDD4hepOK;
31 #endif
32 
33  auto bareBaseName = [&](std::string_view name) {
34  size_t ipos = name.rfind('_');
35  return (isDD4hepOK) ? name.substr(0, ipos) : name;
36  };
37 
38  if (nLevels == kBTLcrystalLevel || isDD4hepOK) {
39  LogDebug("MTDGeom") << bareBaseName(baseNumber.getLevelName(0)) << "[" << baseNumber.getCopyNumber(0) << "], "
40  << bareBaseName(baseNumber.getLevelName(1)) << "[" << baseNumber.getCopyNumber(1) << "], "
41  << bareBaseName(baseNumber.getLevelName(2)) << "[" << baseNumber.getCopyNumber(2) << "], "
42  << bareBaseName(baseNumber.getLevelName(3)) << "[" << baseNumber.getCopyNumber(3) << "], "
43  << bareBaseName(baseNumber.getLevelName(4)) << "[" << baseNumber.getCopyNumber(4) << "], "
44  << bareBaseName(baseNumber.getLevelName(5)) << "[" << baseNumber.getCopyNumber(5) << "], "
45  << bareBaseName(baseNumber.getLevelName(6)) << "[" << baseNumber.getCopyNumber(6) << "], "
46  << bareBaseName(baseNumber.getLevelName(7)) << "[" << baseNumber.getCopyNumber(7) << "], "
47  << bareBaseName(baseNumber.getLevelName(8)) << "[" << baseNumber.getCopyNumber(8) << "]";
48 
49  // barphiflat scenario
50 
51  if (baseNumber.getLevelName(0).find("Timingactive") != std::string_view::npos) {
52  crystal = baseNumber.getCopyNumber(0);
53 
54  modCopy = baseNumber.getCopyNumber(2);
55  rodCopy = baseNumber.getCopyNumber(3);
56 
57  const std::string_view& modName(baseNumber.getLevelName(2)); // name of module volume
58  uint32_t pos = modName.find("Positive");
59 
60  zside = (pos <= modName.size() ? 1 : 0);
61  std::string_view baseName = modName.substr(modName.find(':') + 1);
62 
63  modtyp = ::atoi(&baseName.at(7));
64  if (modtyp == 17) {
65  modtyp = 2;
66  } else if (modtyp == 33) {
67  modtyp = 3;
68  }
69 
70  // error checking
71 
72  if (1 > crystal || 64 < crystal) {
73  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
74  << "****************** Bad crystal number = " << crystal
75  << ", Volume Number = " << baseNumber.getCopyNumber(0);
76  return 0;
77  }
78 
79  if (1 > modtyp || 3 < modtyp) {
80  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
81  << "****************** Bad module name = " << modName
82  << ", Volume Name = " << baseNumber.getLevelName(2);
83  return 0;
84  }
85 
86  if (1 > modCopy || 54 < modCopy) {
87  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
88  << "****************** Bad module copy = " << modCopy
89  << ", Volume Number = " << baseNumber.getCopyNumber(2);
90  return 0;
91  }
92 
93  if (1 > rodCopy || 36 < rodCopy) {
94  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
95  << "****************** Bad rod copy = " << rodCopy
96  << ", Volume Number = " << baseNumber.getCopyNumber(4);
97  return 0;
98  }
99 
100  if (1 < zside) {
101  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
102  << "****************** Bad side = " << zside
103  << ", Volume Name = " << baseNumber.getLevelName(2);
104  return 0;
105  }
106  } else if (baseNumber.getLevelName(0).find("BTLCrystal") != std::string_view::npos) {
107  // v2 or v3 scenario
108 
109  crystal = baseNumber.getCopyNumber(0);
110  modCopy = baseNumber.getCopyNumber(1);
111  runitCopy = baseNumber.getCopyNumber(2);
112  rodCopy = baseNumber.getCopyNumber(3);
113 
114  const std::string_view& rodName(baseNumber.getLevelName(3)); // name of module volume
115  uint32_t pos = rodName.find("Zpos");
116  zside = (pos <= rodName.size() ? 1 : 0);
117 
118  // for negative side swap module numbers betwee sides of the tray, so as to keep the same number for the same phi angle
119  // in the existing model. This introduces a misalignemtn between module number and volume copy for the negative side.
120  if (zside == 0) {
121  modCopy = negModCopy[modCopy - 1];
122  }
123 
124  bool isV2(bareBaseName(baseNumber.getLevelName(0)).back() != 'l');
125 
126 #ifdef EDM_ML_DEBUG
127  LogDebug("MTDGeom") << "BTLNumberingScheme::getUnitID(): isV2 " << isV2;
128 #endif
129 
130  if (isV2) {
131  // V2: the type is embedded in crystal name
132  modtyp = ::atoi(&bareBaseName(baseNumber.getLevelName(2)).back());
133  } else {
134  // V3: build type and RU number per type from global RU number
135  modtyp = globalru2type[runitCopy - 1];
136  runitCopy = globalru2ru[runitCopy - 1];
137  }
138 
139  // error checking
140 
141  if (1 > crystal || BTLDetId::kCrystalsPerModuleV2 < crystal) {
142  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
143  << "****************** Bad crystal number = " << crystal
144  << ", Volume Number = " << baseNumber.getCopyNumber(0);
145  return 0;
146  }
147 
148  if (1 > modtyp || 3 < modtyp) {
149  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
150  << "****************** Bad RU name, Volume Name = "
151  << bareBaseName(baseNumber.getLevelName(2));
152  return 0;
153  }
154 
155  if (1 > modCopy || BTLDetId::kModulesPerRUV2 < modCopy) {
156  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
157  << "****************** Bad module copy = " << modCopy
158  << ", Volume Number = " << baseNumber.getCopyNumber(1);
159  return 0;
160  }
161 
162  if (1 > runitCopy || BTLDetId::kRUPerTypeV2 < runitCopy) {
163  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
164  << "****************** Bad readout unit copy = " << runitCopy
165  << ", Volume Number = " << baseNumber.getCopyNumber(2);
166  return 0;
167  }
168 
169  if (1 > rodCopy || BTLDetId::HALF_ROD < rodCopy) {
170  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
171  << "****************** Bad rod copy = " << rodCopy
172  << ", Volume Number = " << baseNumber.getCopyNumber(3);
173  return 0;
174  }
175 
176  if (1 < zside) {
177  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
178  << "****************** Bad side = " << zside
179  << ", Volume Name = " << baseNumber.getLevelName(3);
180  return 0;
181  }
182  }
183 
184  // all inputs are fine. Go ahead and decode
185 
186  BTLDetId thisBTLdetid(zside, rodCopy, runitCopy, modCopy, modtyp, crystal);
187  intindex = thisBTLdetid.rawId();
188 
189  } else if (nLevels == kBTLmoduleLevel && baseNumber.getLevelName(0).find("BTLModule") != std::string_view::npos) {
190  // v2 scenario, geographicalId per module
191  // for tracking navigation geometry
192  LogDebug("MTDGeom") << bareBaseName(baseNumber.getLevelName(0)) << "[" << baseNumber.getCopyNumber(0) << "], "
193  << bareBaseName(baseNumber.getLevelName(1)) << "[" << baseNumber.getCopyNumber(1) << "], "
194  << bareBaseName(baseNumber.getLevelName(2)) << "[" << baseNumber.getCopyNumber(2) << "], "
195  << bareBaseName(baseNumber.getLevelName(3)) << "[" << baseNumber.getCopyNumber(3) << "], "
196  << bareBaseName(baseNumber.getLevelName(4)) << "[" << baseNumber.getCopyNumber(4) << "], "
197  << bareBaseName(baseNumber.getLevelName(5)) << "[" << baseNumber.getCopyNumber(5) << "], "
198  << bareBaseName(baseNumber.getLevelName(6)) << "[" << baseNumber.getCopyNumber(6) << "], "
199  << bareBaseName(baseNumber.getLevelName(7)) << "[" << baseNumber.getCopyNumber(7) << "]";
200 
201  modCopy = baseNumber.getCopyNumber(0);
202  runitCopy = baseNumber.getCopyNumber(1);
203  rodCopy = baseNumber.getCopyNumber(2);
204 
205  const std::string_view& rodName(baseNumber.getLevelName(2)); // name of module volume
206  uint32_t pos = rodName.find("Zpos");
207  zside = (pos <= rodName.size() ? 1 : 0);
208 
209  // for negative side swap module numbers betwee sides of the tray, so as to keep the same number for the same phi angle
210  // in the existing model. This introduces a misalignemtn between module number and volume copy for the negative side.
211  if (zside == 0) {
212  modCopy = negModCopy[modCopy - 1];
213  }
214 
215  bool isV2(bareBaseName(baseNumber.getLevelName(0)).back() != 'e');
216 
217 #ifdef EDM_ML_DEBUG
218  LogDebug("MTDGeom") << "BTLNumberingScheme::getUnitID(): isV2 " << isV2;
219 #endif
220 
221  if (isV2) {
222  // V2: the type is embedded in crystal name
223  modtyp = ::atoi(&bareBaseName(baseNumber.getLevelName(1)).back());
224  } else {
225  // V3: build type and RU number per type from global RU number
226  modtyp = globalru2type[runitCopy - 1];
227  runitCopy = globalru2ru[runitCopy - 1];
228  }
229 
230  // error checking
231 
232  if (1 > modtyp || 3 < modtyp) {
233  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
234  << "****************** Bad RU name, Volume Name = "
235  << bareBaseName(baseNumber.getLevelName(1));
236  return 0;
237  }
238 
239  if (1 > modCopy || BTLDetId::kModulesPerRUV2 < modCopy) {
240  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
241  << "****************** Bad module copy = " << modCopy
242  << ", Volume Number = " << baseNumber.getCopyNumber(0);
243  return 0;
244  }
245 
246  if (1 > runitCopy || BTLDetId::kRUPerTypeV2 < runitCopy) {
247  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
248  << "****************** Bad readout unit copy = " << runitCopy
249  << ", Volume Number = " << baseNumber.getCopyNumber(1);
250  return 0;
251  }
252 
253  if (1 > rodCopy || BTLDetId::HALF_ROD < rodCopy) {
254  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
255  << "****************** Bad rod copy = " << rodCopy
256  << ", Volume Number = " << baseNumber.getCopyNumber(2);
257  return 0;
258  }
259 
260  if (1 < zside) {
261  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
262  << "****************** Bad side = " << zside
263  << ", Volume Name = " << baseNumber.getLevelName(2);
264  return 0;
265  }
266 
267  // all inputs are fine. Go ahead and decode
268 
269  BTLDetId thisBTLdetid(zside, rodCopy, runitCopy, modCopy, modtyp, 0);
270  intindex = thisBTLdetid.geographicalId(BTLDetId::CrysLayout::v2).rawId();
271 
272  } else {
273  edm::LogWarning("MTDGeom") << "BTLNumberingScheme::getUnitID(): "
274  << "Not enough levels found in MTDBaseNumber ( " << nLevels
275  << ") or not correct path. Returning 0";
276  return 0;
277  }
278 
279  LogDebug("MTDGeom") << "BTL Numbering scheme: "
280  << " zside = " << zside << " rod = " << rodCopy << " modtyp = " << modtyp << " RU = " << runitCopy
281  << " module = " << modCopy << " crystal = " << crystal << " Raw Id = " << intindex << "\n"
282  << BTLDetId(intindex);
283 
284  return intindex;
285 }
int getCopyNumber(int level) const
static constexpr std::array< uint32_t, BTLDetId::kRUPerTypeV2 *BTLDetId::kCrystalTypes > globalru2ru
static constexpr std::array< uint32_t, BTLDetId::kModulesPerRUV2 > negModCopy
int getLevels() const
int zside(DetId const &)
std::string_view const & getLevelName(int level) const
static constexpr uint32_t HALF_ROD
range constants, need two sets for the time being (one for tiles and one for bars) ...
Definition: BTLDetId.h:31
static constexpr uint32_t kModulesPerRUV2
Definition: BTLDetId.h:35
static constexpr uint32_t kBTLcrystalLevel
static constexpr uint32_t kRUPerTypeV2
Definition: BTLDetId.h:34
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:19
static constexpr uint32_t kCrystalsPerModuleV2
Definition: BTLDetId.h:36
static constexpr uint32_t kBTLmoduleLevel
Log< level::Warning, false > LogWarning
static constexpr std::array< uint32_t, BTLDetId::kRUPerTypeV2 *BTLDetId::kCrystalTypes > globalru2type
#define LogDebug(id)

Member Data Documentation

◆ globalru2ru

constexpr std::array<uint32_t, BTLDetId::kRUPerTypeV2 * BTLDetId::kCrystalTypes> BTLNumberingScheme::globalru2ru
static
Initial value:
{
{1, 2, 1, 2, 1, 2}}

Definition at line 18 of file BTLNumberingScheme.h.

Referenced by getUnitID().

◆ globalru2type

constexpr std::array<uint32_t, BTLDetId::kRUPerTypeV2 * BTLDetId::kCrystalTypes> BTLNumberingScheme::globalru2type
static
Initial value:
{
{1, 1, 2, 2, 3, 3}}

Definition at line 16 of file BTLNumberingScheme.h.

Referenced by getUnitID().

◆ kBTLcrystalLevel

constexpr uint32_t BTLNumberingScheme::kBTLcrystalLevel = 9
static

Definition at line 9 of file BTLNumberingScheme.h.

Referenced by getUnitID().

◆ kBTLmoduleLevel

constexpr uint32_t BTLNumberingScheme::kBTLmoduleLevel = 8
static

Definition at line 10 of file BTLNumberingScheme.h.

Referenced by getUnitID().

◆ negModCopy

constexpr std::array<uint32_t, BTLDetId::kModulesPerRUV2> BTLNumberingScheme::negModCopy
static
Initial value:
{
{3, 2, 1, 6, 5, 4, 9, 8, 7, 12, 11, 10, 15, 14, 13, 18, 17, 16, 21, 20, 19, 24, 23, 22}}

Definition at line 12 of file BTLNumberingScheme.h.

Referenced by getUnitID().