File indexing completed on 2024-05-10 02:21:19
0001
0002
0003
0004
0005 #include "SimG4CMS/Forward/interface/CastorNumberingScheme.h"
0006 #include "DataFormats/HcalDetId/interface/HcalCastorDetId.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008
0009 #include <CLHEP/Units/SystemOfUnits.h>
0010 #include "G4LogicalVolumeStore.hh"
0011 #include <iostream>
0012
0013
0014
0015 CastorNumberingScheme::CastorNumberingScheme()
0016 : lvCASTFar(nullptr),
0017 lvCASTNear(nullptr),
0018 lvCAST(nullptr),
0019 lvCAES(nullptr),
0020 lvCEDS(nullptr),
0021 lvCAHS(nullptr),
0022 lvCHDS(nullptr),
0023 lvCAER(nullptr),
0024 lvCEDR(nullptr),
0025 lvCAHR(nullptr),
0026 lvCHDR(nullptr),
0027 lvC3EF(nullptr),
0028 lvC3HF(nullptr),
0029 lvC4EF(nullptr),
0030 lvC4HF(nullptr) {
0031 const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
0032 std::vector<lvp>::const_iterator lvcite;
0033 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
0034 if (strcmp(((*lvcite)->GetName()).c_str(), "CASTFar") == 0)
0035 lvCASTFar = (*lvcite);
0036 if (strcmp(((*lvcite)->GetName()).c_str(), "CASTNear") == 0)
0037 lvCASTNear = (*lvcite);
0038 if (strcmp(((*lvcite)->GetName()).c_str(), "CAST") == 0)
0039 lvCAST = (*lvcite);
0040 if (strcmp(((*lvcite)->GetName()).c_str(), "CAES") == 0)
0041 lvCAES = (*lvcite);
0042 if (strcmp(((*lvcite)->GetName()).c_str(), "CEDS") == 0)
0043 lvCEDS = (*lvcite);
0044 if (strcmp(((*lvcite)->GetName()).c_str(), "CAHS") == 0)
0045 lvCAHS = (*lvcite);
0046 if (strcmp(((*lvcite)->GetName()).c_str(), "CHDS") == 0)
0047 lvCHDS = (*lvcite);
0048 if (strcmp(((*lvcite)->GetName()).c_str(), "CAER") == 0)
0049 lvCAER = (*lvcite);
0050 if (strcmp(((*lvcite)->GetName()).c_str(), "CEDR") == 0)
0051 lvCEDR = (*lvcite);
0052 if (strcmp(((*lvcite)->GetName()).c_str(), "CAHR") == 0)
0053 lvCAHR = (*lvcite);
0054 if (strcmp(((*lvcite)->GetName()).c_str(), "CHDR") == 0)
0055 lvCHDR = (*lvcite);
0056 if (strcmp(((*lvcite)->GetName()).c_str(), "C3EF") == 0)
0057 lvC3EF = (*lvcite);
0058 if (strcmp(((*lvcite)->GetName()).c_str(), "C3HF") == 0)
0059 lvC3HF = (*lvcite);
0060 if (strcmp(((*lvcite)->GetName()).c_str(), "C4EF") == 0)
0061 lvC4EF = (*lvcite);
0062 if (strcmp(((*lvcite)->GetName()).c_str(), "C4HF") == 0)
0063 lvC4HF = (*lvcite);
0064 }
0065 #ifdef EDM_ML_DEBUG
0066 edm::LogVerbatim("ForwardSim") << "Creating CastorNumberingScheme";
0067 LogDebug("ForwardSim") << "CastorNumberingScheme:: LogicalVolume pointers\n"
0068 << lvCASTFar << " for CASTFar; " << lvCASTNear << " for CASTNear; " << lvCAST << " for CAST; "
0069 << lvCAES << " for CAES; " << lvCEDS << " for CEDS; " << lvCAHS << " for CAHS; " << lvCHDS
0070 << " for CHDS; " << lvCAER << " for CAER; " << lvCEDR << " for CEDR; " << lvCAHR
0071 << " for CAHR; " << lvCHDR << " for CHDR; " << lvC3EF << " for C3EF; " << lvC3HF
0072 << " for C3HF; " << lvC4EF << " for C4EF; " << lvC4HF << " for C4HF.";
0073
0074 LogDebug("ForwardSim") << "Call to init CastorNumberingScheme\n";
0075 for (int mod = 0; mod < 15; mod++)
0076 for (int sec = 0; sec < 17; sec++) {
0077 HcalCastorDetId castorId = HcalCastorDetId(false, sec, mod);
0078 LogDebug("ForwardSim") << "Mod: " << mod << " Sec: " << sec << " Id: " << castorId.rawId() << "\n";
0079 }
0080
0081 #endif
0082 }
0083
0084 CastorNumberingScheme::~CastorNumberingScheme() {}
0085
0086 uint32_t CastorNumberingScheme::getUnitID(const G4Step* aStep) const {
0087 uint32_t index = 0;
0088 int level, copyno[20];
0089 lvp lvs[20];
0090 detectorLevel(aStep, level, copyno, lvs);
0091
0092 #ifdef EDM_ML_DEBUG
0093 LogDebug("ForwardSim") << "CastorNumberingScheme number of levels= " << level;
0094 #endif
0095
0096 if (level > 0) {
0097 int zside = 0;
0098 int sector = 0;
0099 int module = 0;
0100
0101 bool farSide = false;
0102 int castorGeoVersion = 0;
0103
0104
0105 for (int ich = 0; ich < level; ich++) {
0106 if (lvs[ich] == lvCAST) {
0107
0108 assert(1 <= copyno[ich] && copyno[ich] <= 3);
0109 zside = copyno[ich] == 1 ? 1 : 2;
0110 }
0111 else if (lvs[ich] == lvCASTFar || lvs[ich] == lvCASTNear) {
0112 castorGeoVersion = 1;
0113 if (lvs[ich] == lvCASTFar)
0114 farSide = true;
0115 } else if (lvs[ich] == lvCAES || lvs[ich] == lvCEDS || lvs[ich] == lvCAHS || lvs[ich] == lvCHDS) {
0116
0117 int copyn = copyno[ich];
0118 if (castorGeoVersion == 1) {
0119
0120
0121 if (farSide) {
0122 if (copyn < 3)
0123 copyn += 6;
0124 else
0125 copyn -= 2;
0126 } else {
0127 copyn += 2;
0128 }
0129 }
0130 if (copyn < 5)
0131 sector = 5 - copyn;
0132 else
0133 sector = 13 - copyn;
0134 } else if (lvs[ich] == lvCAER || lvs[ich] == lvCEDR) {
0135
0136 module = copyno[ich];
0137 } else if (lvs[ich] == lvCAHR || lvs[ich] == lvCHDR) {
0138
0139 module = copyno[ich] + 2;
0140 } else if (lvs[ich] == lvC3EF || lvs[ich] == lvC3HF) {
0141
0142 sector = sector * 2;
0143 } else if (lvs[ich] == lvC4EF || lvs[ich] == lvC4HF) {
0144
0145 sector = sector * 2 - 1;
0146 }
0147
0148 #ifdef EDM_ML_DEBUG
0149 LogDebug("ForwardSim") << "CastorNumberingScheme "
0150 << "ich = " << ich << "copyno = " << copyno[ich] << "name = " << lvs[ich]->GetName();
0151 #endif
0152 }
0153
0154
0155
0156
0157
0158
0159
0160 bool true_for_positive_eta = false;
0161 if (zside == 1)
0162 true_for_positive_eta = true;
0163
0164 HcalCastorDetId castorId = HcalCastorDetId(true_for_positive_eta, sector, module);
0165 index = castorId.rawId();
0166
0167 #ifdef EDM_ML_DEBUG
0168 uint32_t intindex = 0;
0169 intindex = packIndex(zside, sector, module);
0170 LogDebug("ForwardSim") << "CastorNumberingScheme: "
0171 << " zside " << zside << " module " << module << " sector " << sector << " UnitID 0x"
0172 << std::hex << intindex << std::dec << " index: " << index;
0173 #endif
0174 }
0175 return index;
0176 }
0177
0178 uint32_t CastorNumberingScheme::packIndex(int z, int sector, int module) {
0179
0180
0181
0182
0183
0184
0185
0186
0187 uint32_t idx = ((z - 1) & 1) << 8;
0188 idx += (sector & 15) << 4;
0189 idx += (module & 15);
0190 return idx;
0191 }
0192
0193 void CastorNumberingScheme::unpackIndex(const uint32_t& idx, int& z, int& sector, int& module) {
0194
0195
0196
0197
0198
0199
0200
0201 z = (idx >> 8) & 1;
0202 z += 1;
0203 sector = (idx >> 4) & 15;
0204 module = (idx & 15);
0205 }
0206
0207 void CastorNumberingScheme::detectorLevel(const G4Step* aStep, int& level, int* copyno, lvp* lvs) const {
0208
0209 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
0210 level = 0;
0211 if (touch)
0212 level = ((touch->GetHistoryDepth()) + 1);
0213 if (level > 0) {
0214 for (int ii = 0; ii < level; ii++) {
0215 int i = level - ii - 1;
0216 lvs[ii] = touch->GetVolume(i)->GetLogicalVolume();
0217 copyno[ii] = touch->GetReplicaNumber(i);
0218 }
0219 }
0220 }