Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-25 23:57:06

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "Geometry/CaloGeometry/interface/CaloGenericDetId.h"
0003 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
0004 #include "Geometry/ForwardGeometry/interface/ZdcGeometry.h"
0005 #include "Geometry/ForwardGeometry/interface/IdealZDCTrapezoid.h"
0006 #include "Geometry/ForwardGeometry/interface/ZdcHardcodeGeometryData.h"
0007 #include <algorithm>
0008 
0009 #include <Math/Transform3D.h>
0010 #include <Math/EulerAngles.h>
0011 
0012 #include <algorithm>
0013 
0014 typedef CaloCellGeometry::Pt3D Pt3D;
0015 typedef CaloCellGeometry::Pt3DVec Pt3DVec;
0016 typedef CaloCellGeometry::Tr3D Tr3D;
0017 
0018 typedef CaloSubdetectorGeometry::CCGFloat CCGFloat;
0019 
0020 //#define EDM_ML_DEBUG
0021 
0022 ZdcGeometry::ZdcGeometry()
0023     : k_NumberOfCellsForCornersN(HcalZDCDetId::kSizeForDenseIndexingRun1),
0024       theTopology(nullptr),
0025       lastReqDet_(DetId::Detector(0)),
0026       lastReqSubdet_(0),
0027       m_ownsTopology(true),
0028       m_cellVec(k_NumberOfCellsForCornersN) {
0029   edm::LogWarning("HCalGeom") << "ZdcGeometry::Wrong constructor called";
0030 }
0031 
0032 ZdcGeometry::ZdcGeometry(const ZdcTopology* topology)
0033     : k_NumberOfCellsForCornersN(topology->kSizeForDenseIndexing()),
0034       theTopology(topology),
0035       lastReqDet_(DetId::Detector(0)),
0036       lastReqSubdet_(0),
0037       m_ownsTopology(false),
0038       m_cellVec(k_NumberOfCellsForCornersN) {}
0039 
0040 ZdcGeometry::~ZdcGeometry() {
0041   if (m_ownsTopology)
0042     delete theTopology;
0043 }
0044 /*
0045 DetId ZdcGeometry::getClosestCell(const GlobalPoint& r) const
0046 {
0047    DetId returnId ( 0 ) ;
0048    const std::vector<DetId>& detIds ( getValidDetIds() ) ;
0049    for( std::vector<DetId>::const_iterator it ( detIds.begin() ) ;
0050     it != detIds.end(); ++it )
0051    {
0052       auto cell = ( getGeometry( *it ) ) ;
0053       if( 0 != cell &&
0054       cell->inside( r ) )
0055       {
0056      returnId = *it ;
0057      break ;
0058       }
0059    }
0060    return returnId ;
0061 }
0062 */
0063 unsigned int ZdcGeometry::alignmentTransformIndexLocal(const DetId& id) {
0064   const CaloGenericDetId gid(id);
0065   assert(gid.isZDC());
0066 
0067   return (0 > HcalZDCDetId(id).zside() ? 0 : 1);
0068 }
0069 
0070 unsigned int ZdcGeometry::alignmentTransformIndexGlobal(const DetId& /*id*/) { return (unsigned int)DetId::Calo - 1; }
0071 
0072 void ZdcGeometry::localCorners(Pt3DVec& lc, const CCGFloat* pv, unsigned int /*i*/, Pt3D& ref) {
0073   IdealZDCTrapezoid::localCorners(lc, pv, ref);
0074 }
0075 
0076 void ZdcGeometry::newCell(const GlobalPoint& f1,
0077                           const GlobalPoint& /*f2*/,
0078                           const GlobalPoint& /*f3*/,
0079                           const CCGFloat* parm,
0080                           const DetId& detId) {
0081   const CaloGenericDetId cgid(detId);
0082 #ifdef EDM_ML_DEBUG
0083   edm::LogVerbatim("ZDCGeom") << "ZDCGeometry " << HcalZDCDetId(detId) << " Generic ID " << std::hex << cgid.rawId()
0084                               << std::dec << " ZDC? " << cgid.isZDC();
0085 #endif
0086   assert(cgid.isZDC());
0087 
0088   const unsigned int di(theTopology->detId2DenseIndex(detId));
0089 
0090   m_cellVec[di] = IdealZDCTrapezoid(f1, cornersMgr(), parm);
0091   addValidID(detId);
0092 }
0093 
0094 CaloCellGeometryPtr ZdcGeometry::getGeometryRawPtr(uint32_t index) const {
0095   // Modify the RawPtr class
0096   if (m_cellVec.size() <= index)
0097     return CaloCellGeometryPtr();
0098   const CaloCellGeometry* cell(&m_cellVec[index]);
0099   return CaloCellGeometryPtr(((cell == nullptr) || (nullptr == cell->param())) ? nullptr : cell);
0100 }
0101 
0102 void ZdcGeometry::getSummary(CaloSubdetectorGeometry::TrVec& tVec,
0103                              CaloSubdetectorGeometry::IVec& iVec,
0104                              CaloSubdetectorGeometry::DimVec& dVec,
0105                              CaloSubdetectorGeometry::IVec& /*dins*/) const {
0106   tVec.reserve(m_validIds.size() * numberOfTransformParms());
0107   iVec.reserve(numberOfShapes() == 1 ? 1 : m_validIds.size());
0108   dVec.reserve(numberOfShapes() * numberOfParametersPerShape());
0109 
0110   for (const auto& pv : parVecVec()) {
0111     for (float iv : pv) {
0112       dVec.emplace_back(iv);
0113     }
0114   }
0115 
0116   for (uint32_t i(0); i != m_validIds.size(); ++i) {
0117     Tr3D tr;
0118     auto ptr = cellGeomPtr(i);
0119 #ifdef EDM_ML_DEBUG
0120     edm::LogVerbatim("ZDCGeom") << "ZDCGeometry:Summary " << i << ":" << HcalZDCDetId::kSizeForDenseIndexingRun1
0121                                 << " Pointer " << ptr << ":" << (nullptr != ptr);
0122 #endif
0123     if (i < static_cast<int32_t>(HcalZDCDetId::kSizeForDenseIndexingRun1))
0124       assert(nullptr != ptr);
0125     if (ptr != nullptr) {
0126       ptr->getTransform(tr, (Pt3DVec*)nullptr);
0127 
0128       if (Tr3D() == tr) {  // for preshower there is no rotation
0129         const GlobalPoint& gp(ptr->getPosition());
0130         tr = HepGeom::Translate3D(gp.x(), gp.y(), gp.z());
0131       }
0132 
0133       const CLHEP::Hep3Vector tt(tr.getTranslation());
0134       tVec.emplace_back(tt.x());
0135       tVec.emplace_back(tt.y());
0136       tVec.emplace_back(tt.z());
0137       if (6 == numberOfTransformParms()) {
0138         const CLHEP::HepRotation rr(tr.getRotation());
0139         const ROOT::Math::Transform3D rtr(
0140             rr.xx(), rr.xy(), rr.xz(), tt.x(), rr.yx(), rr.yy(), rr.yz(), tt.y(), rr.zx(), rr.zy(), rr.zz(), tt.z());
0141         ROOT::Math::EulerAngles ea;
0142         rtr.GetRotation(ea);
0143         tVec.emplace_back(ea.Phi());
0144         tVec.emplace_back(ea.Theta());
0145         tVec.emplace_back(ea.Psi());
0146       }
0147 
0148       const CCGFloat* par(ptr->param());
0149 
0150       unsigned int ishape(9999);
0151       for (unsigned int ivv(0); ivv != parVecVec().size(); ++ivv) {
0152         bool ok(true);
0153         const CCGFloat* pv(&(*parVecVec()[ivv].begin()));
0154         for (unsigned int k(0); k != numberOfParametersPerShape(); ++k) {
0155           ok = ok && (fabs(par[k] - pv[k]) < 1.e-6);
0156         }
0157         if (ok) {
0158           ishape = ivv;
0159           break;
0160         }
0161       }
0162       assert(9999 != ishape);
0163 
0164       const unsigned int nn((numberOfShapes() == 1) ? (unsigned int)1 : m_validIds.size());
0165       if (iVec.size() < nn)
0166         iVec.emplace_back(ishape);
0167     }
0168   }
0169 }