Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-21 01:40:02

0001 #include "RecoLocalTracker/Phase2TrackerRecHits/interface/Phase2StripCPE.h"
0002 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 
0005 Phase2StripCPE::Phase2StripCPE(edm::ParameterSet& conf,
0006                                const MagneticField& magf,
0007                                const TrackerGeometry& geom,
0008                                const SiPhase2OuterTrackerLorentzAngle& LorentzAngle)
0009     : magfield_(magf),
0010       geom_(geom),
0011       lorentzAngleMap_(LorentzAngle),
0012       tanLorentzAnglePerTesla_(conf.getParameter<double>("TanLorentzAnglePerTesla")) {
0013   use_LorentzAngle_DB_ = conf.getParameter<bool>("LorentzAngle_DB");
0014   fillParam();
0015 }
0016 
0017 void Phase2StripCPE::fillPSetDescription(edm::ParameterSetDescription& desc) {
0018   desc.add<double>("TanLorentzAnglePerTesla", 0.07);
0019   desc.add<bool>("LorentzAngle_DB", true);
0020 }
0021 
0022 Phase2StripCPE::LocalValues Phase2StripCPE::localParameters(const Phase2TrackerCluster1D& cluster,
0023                                                             const GeomDetUnit& detunit) const {
0024   auto const& p = m_Params[detunit.index() - m_off];
0025   auto const& topo = *p.topology;
0026   float ix = cluster.center() - 0.5f * p.coveredStrips;
0027   float iy = float(cluster.column()) + 0.5f;  // halfway the column
0028 
0029   LocalPoint lp(topo.localX(ix), topo.localY(iy), 0);  // x, y, z
0030 
0031   return std::make_pair(lp, p.localErr);
0032 }
0033 
0034 LocalVector Phase2StripCPE::driftDirection(const Phase2TrackerGeomDetUnit& det) const {
0035   LocalVector lbfield = (det.surface()).toLocal(magfield_.inTesla(det.surface().position()));
0036 
0037   float langle =
0038       use_LorentzAngle_DB_ ? lorentzAngleMap_.getLorentzAngle(det.geographicalId().rawId()) : tanLorentzAnglePerTesla_;
0039 
0040   float dir_x = -langle * lbfield.y();
0041   float dir_y = langle * lbfield.x();
0042   float dir_z = 1.f;  // E field always in z direction
0043 
0044   return LocalVector(dir_x, dir_y, dir_z);
0045 }
0046 
0047 void Phase2StripCPE::fillParam() {
0048   // in phase 2 they are all pixel topologies...
0049   auto const& dus = geom_.detUnits();
0050   m_off = dus.size();
0051   // skip Barrel and Foward pixels...
0052   for (unsigned int i = 3; i < 7; ++i) {
0053     LogDebug("LookingForFirstPhase2OT") << " Subdetector " << i << " GeomDetEnumerator "
0054                                         << GeomDetEnumerators::tkDetEnum[i] << " offset "
0055                                         << geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i]) << std::endl;
0056     if (geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i]) != dus.size()) {
0057       if (geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i]) < m_off)
0058         m_off = geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i]);
0059     }
0060   }
0061   LogDebug("LookingForFirstPhase2OT") << " Chosen offset: " << m_off;
0062 
0063   m_Params.resize(dus.size() - m_off);
0064   // very very minimal, for sure it will need to expand...
0065   for (auto i = m_off; i != dus.size(); ++i) {
0066     auto& p = m_Params[i - m_off];
0067 
0068     const Phase2TrackerGeomDetUnit& det = (const Phase2TrackerGeomDetUnit&)(*dus[i]);
0069     assert(det.index() == int(i));
0070     p.topology = &det.specificTopology();
0071 
0072     auto pitch_x = p.topology->pitch().first;
0073     auto pitch_y = p.topology->pitch().second;
0074 
0075     // see https://github.com/cms-sw/cmssw/blob/CMSSW_8_1_X/RecoLocalTracker/SiStripRecHitConverter/src/StripCPE.cc
0076     auto thickness = det.specificSurface().bounds().thickness();
0077     auto drift = driftDirection(det) * thickness;
0078     auto lvec = drift + LocalVector(0, 0, -thickness);
0079     p.coveredStrips = lvec.x() / pitch_x;  // simplifies wrt Phase0 tracker because only rectangular modules
0080 
0081     constexpr float o12 = 1. / 12;
0082     p.localErr = LocalError(o12 * pitch_x * pitch_x, 0, o12 * pitch_y * pitch_y);  // e2_xx, e2_xy, e2_yy
0083   }
0084 }