Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:14:47

0001 #include "Geometry/HcalCommonData/interface/HcalHitRelabeller.h"
0002 #include "DataFormats/HcalDetId/interface/HcalTestNumbering.h"
0003 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0004 
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 
0007 //#define EDM_ML_DEBUG
0008 
0009 HcalHitRelabeller::HcalHitRelabeller(bool nd) : theRecNumber(nullptr), neutralDensity_(nd) {
0010 #ifdef EDM_ML_DEBUG
0011   edm::LogVerbatim("HcalSim") << "HcalHitRelabeller initialized with"
0012                               << " neutralDensity " << neutralDensity_;
0013 #endif
0014 }
0015 
0016 void HcalHitRelabeller::process(std::vector<PCaloHit>& hcalHits) {
0017   if (theRecNumber) {
0018 #ifdef EDM_ML_DEBUG
0019     int ii(0);
0020 #endif
0021     for (auto& hcalHit : hcalHits) {
0022 #ifdef EDM_ML_DEBUG
0023       edm::LogVerbatim("HcalSim") << "Hit[" << ii << "] " << std::hex << hcalHit.id() << std::dec << " Neutral density "
0024                                   << neutralDensity_;
0025 #endif
0026       uint32_t hid;
0027       int det, z, depth, eta, phi, layer;
0028       HcalTestNumbering::unpackHcalIndex(hcalHit.id(), det, z, depth, eta, phi, layer);
0029       if ((det == 2) && (layer == 2) && (eta == 18))
0030         depth = 2;
0031       hid = HcalTestNumbering::packHcalIndex(det, z, depth, eta, phi, layer);
0032       double wt = (neutralDensity_) ? (energyWt(hid)) : 1.0;
0033       double energy = wt * (hcalHit.energy());
0034       hcalHit.setEnergy(energy);
0035       uint32_t newid = relabel(hid).rawId();
0036 #ifdef EDM_ML_DEBUG
0037       edm::LogVerbatim("HcalSim") << "Hit " << ii << " out of " << hcalHits.size() << " " << std::hex << newid
0038                                   << std::dec << " E " << energy << " wt " << wt;
0039 #endif
0040       hcalHit.setID(newid);
0041 #ifdef EDM_ML_DEBUG
0042       edm::LogVerbatim("HcalSim") << "Final Setting::subdet: " << HcalDetId(newid).subdet()
0043                                   << " z: " << HcalDetId(newid).zside() << " depth: " << HcalDetId(newid).depth()
0044                                   << " ieta: " << HcalDetId(newid).ietaAbs() << " iphi: " << HcalDetId(newid).iphi()
0045                                   << " wt " << wt;
0046       ++ii;
0047 #endif
0048     }
0049   } else {
0050     edm::LogWarning("HcalSim") << "HcalHitRelabeller: no valid HcalDDDRecConstants";
0051   }
0052 }
0053 
0054 void HcalHitRelabeller::setGeometry(const HcalDDDRecConstants*& recNum) { theRecNumber = recNum; }
0055 
0056 DetId HcalHitRelabeller::relabel(const uint32_t testId) const {
0057   return HcalHitRelabeller::relabel(testId, theRecNumber);
0058 }
0059 
0060 DetId HcalHitRelabeller::relabel(const uint32_t testId, const HcalDDDRecConstants* theRecNumber) {
0061 #ifdef EDM_ML_DEBUG
0062   edm::LogVerbatim("HcalSim") << "Enter HcalHitRelabeller::relabel";
0063 #endif
0064   HcalDetId hid;
0065   int det, z, depth, eta, phi, layer, sign;
0066   HcalTestNumbering::unpackHcalIndex(testId, det, z, depth, eta, phi, layer);
0067 
0068   sign = (z == 0) ? (-1) : (1);
0069   HcalDDDRecConstants::HcalID id = theRecNumber->getHCID(det, sign * eta, phi, layer, depth);
0070   int depth0 = ((det == 1) && (layer == 2) && (eta == 15)) ? depth : id.depth;
0071 
0072   if (id.subdet == int(HcalBarrel)) {
0073     hid = HcalDetId(HcalBarrel, sign * id.eta, id.phi, depth0);
0074   } else if (id.subdet == int(HcalEndcap)) {
0075     hid = HcalDetId(HcalEndcap, sign * id.eta, id.phi, depth0);
0076   } else if (id.subdet == int(HcalOuter)) {
0077     hid = HcalDetId(HcalOuter, sign * id.eta, id.phi, depth0);
0078   } else if (id.subdet == int(HcalForward)) {
0079     hid = HcalDetId(HcalForward, sign * id.eta, id.phi, depth0);
0080   }
0081 #ifdef EDM_ML_DEBUG
0082   edm::LogVerbatim("HcalSim") << "Initial indices:"
0083                               << "subdet: " << det << " "
0084                               << "z: " << z << " "
0085                               << "depth: " << depth << " "
0086                               << "ieta: " << eta << " "
0087                               << "iphi: " << phi << " "
0088                               << "layer: " << layer << " new HcalDetId -> hex.RawID = " << std::hex << hid.rawId()
0089                               << std::dec << " subdet, z, depth, eta, phi = " << det << " " << z << " " << id.depth
0090                               << " " << id.eta << " " << id.phi << " ---> " << hid;
0091 #endif
0092   return hid;
0093 }
0094 
0095 double HcalHitRelabeller::energyWt(const uint32_t testId) const {
0096   int det, z, depth, eta, phi, layer;
0097   HcalTestNumbering::unpackHcalIndex(testId, det, z, depth, eta, phi, layer);
0098   int zside = (z == 0) ? (-1) : (1);
0099   double wt = ((((det == 1) && (layer <= 1)) || ((det == 2) && (layer <= 2))) && (depth == 1))
0100                   ? theRecNumber->getLayer0Wt(det, phi, zside)
0101                   : 1.0;
0102 #ifdef EDM_ML_DEBUG
0103   edm::LogVerbatim("HcalSim") << "EnergyWT::det: " << det << " z: " << z << ":" << zside << " depth: " << depth
0104                               << " ieta: " << eta << " iphi: " << phi << " layer: " << layer << " wt " << wt;
0105 #endif
0106   return wt;
0107 }