Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:35

0001 #include <memory>
0002 #include <string>
0003 #include <iostream>
0004 #include <TMath.h>
0005 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTLA.h"
0006 
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/ConsumesCollector.h"
0009 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0010 #include "DataFormats/GeometryVector/interface/LocalVector.h"
0011 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0012 #include "Geometry/CommonDetUnit/interface/GluedGeomDet.h"
0013 #include "Geometry/CommonTopologies/interface/Topology.h"
0014 #include "Geometry/CommonTopologies/interface/StripTopology.h"
0015 #include "DataFormats/TrackReco/interface/Track.h"
0016 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0017 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0018 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
0019 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0020 #include "TrackingTools/PatternTools/interface/Trajectory.h"
0021 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0022 
0023 using namespace std;
0024 SiStripFineDelayTLA::SiStripFineDelayTLA(edm::ParameterSet const& conf, edm::ConsumesCollector iC)
0025     : conf_(conf), tkGeomToken_(iC.esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()) {}
0026 
0027 void SiStripFineDelayTLA::init(const edm::Event& e, const edm::EventSetup& es) { tracker = &es.getData(tkGeomToken_); }
0028 
0029 // Virtual destructor needed.
0030 SiStripFineDelayTLA::~SiStripFineDelayTLA() {}
0031 
0032 std::vector<std::pair<std::pair<DetId, LocalPoint>, float> > SiStripFineDelayTLA::findtrackangle(
0033     const std::vector<Trajectory>& trajVec) {
0034   if (!trajVec.empty()) {
0035     return findtrackangle(trajVec.front());
0036   }
0037   std::vector<std::pair<std::pair<DetId, LocalPoint>, float> > hitangleassociation;
0038   return hitangleassociation;
0039 }
0040 
0041 std::vector<std::pair<std::pair<DetId, LocalPoint>, float> > SiStripFineDelayTLA::findtrackangle(
0042     const Trajectory& traj) {
0043   std::vector<std::pair<std::pair<DetId, LocalPoint>, float> > hitangleassociation;
0044   std::vector<TrajectoryMeasurement> TMeas = traj.measurements();
0045   std::vector<TrajectoryMeasurement>::iterator itm;
0046   for (itm = TMeas.begin(); itm != TMeas.end(); itm++) {
0047     TrajectoryStateOnSurface tsos = itm->updatedState();
0048     auto thit = itm->recHit();
0049     const SiStripMatchedRecHit2D* matchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>((*thit).hit());
0050     const SiStripRecHit2D* hit = dynamic_cast<const SiStripRecHit2D*>((*thit).hit());
0051     LocalVector trackdirection = tsos.localDirection();
0052     if (matchedhit) {  //if matched hit...
0053       const GluedGeomDet* gdet = static_cast<const GluedGeomDet*>(tracker->idToDet(matchedhit->geographicalId()));
0054       GlobalVector gtrkdir = gdet->toGlobal(trackdirection);
0055       // trackdirection on mono det
0056       // this the pointer to the mono hit of a matched hit
0057       const SiStripRecHit2D monohit = matchedhit->monoHit();
0058       const GeomDetUnit* monodet = gdet->monoDet();
0059       LocalVector monotkdir = monodet->toLocal(gtrkdir);
0060       if (monotkdir.z() != 0) {
0061         // the local angle (mono)
0062         float localpitch = static_cast<const StripTopology&>(monodet->topology()).localPitch(tsos.localPosition());
0063         float thickness = ((((((monohit.geographicalId()) >> 25) & 0x7f) == 0xd) ||
0064                             ((((monohit.geographicalId()) >> 25) & 0x7f) == 0xe)) &&
0065                            ((((monohit.geographicalId()) >> 5) & 0x7) > 4))
0066                               ? 0.0500
0067                               : 0.0320;
0068         float angle = computeAngleCorr(monotkdir, localpitch, thickness);
0069         hitangleassociation.push_back(make_pair(make_pair(monohit.geographicalId(), monohit.localPosition()), angle));
0070         // trackdirection on stereo det
0071         // this the pointer to the stereo hit of a matched hit
0072         const SiStripRecHit2D stereohit = matchedhit->stereoHit();
0073         const GeomDetUnit* stereodet = gdet->stereoDet();
0074         LocalVector stereotkdir = stereodet->toLocal(gtrkdir);
0075         if (stereotkdir.z() != 0) {
0076           // the local angle (stereo)
0077           float localpitch = static_cast<const StripTopology&>(stereodet->topology()).localPitch(tsos.localPosition());
0078           float thickness = ((((((stereohit.geographicalId()) >> 25) & 0x7f) == 0xd) ||
0079                               ((((stereohit.geographicalId()) >> 25) & 0x7f) == 0xe)) &&
0080                              ((((stereohit.geographicalId()) >> 5) & 0x7) > 4))
0081                                 ? 0.0500
0082                                 : 0.0320;
0083           float angle = computeAngleCorr(stereotkdir, localpitch, thickness);
0084           hitangleassociation.push_back(
0085               make_pair(make_pair(stereohit.geographicalId(), stereohit.localPosition()), angle));
0086         }
0087       }
0088     } else if (hit) {
0089       const GeomDetUnit* det = tracker->idToDet(hit->geographicalId());
0090       //  hit= pointer to the rechit
0091       if (trackdirection.z() != 0) {
0092         // the local angle (single hit)
0093         float localpitch = static_cast<const StripTopology&>(det->topology()).localPitch(tsos.localPosition());
0094         float thickness =
0095             ((((((hit->geographicalId()) >> 25) & 0x7f) == 0xd) || ((((hit->geographicalId()) >> 25) & 0x7f) == 0xe)) &&
0096              ((((hit->geographicalId()) >> 5) & 0x7) > 4))
0097                 ? 0.0500
0098                 : 0.0320;
0099         float angle = computeAngleCorr(trackdirection, localpitch, thickness);
0100         hitangleassociation.push_back(make_pair(make_pair(hit->geographicalId(), hit->localPosition()), angle));
0101       }
0102     }
0103   }
0104   return hitangleassociation;
0105 }
0106 
0107 double SiStripFineDelayTLA::computeAngleCorr(const LocalVector& v, double pitch, double thickness) {
0108   double v_xy = sqrt(v.x() * v.x() + v.y() * v.y());
0109   double L = fabs(thickness * v_xy / v.z());
0110   double Lmax = fabs(pitch / v.x() * v_xy);
0111   if (L < Lmax) {
0112     LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax << " Signal contained in strip. Correction is "
0113                                      << v.z() / v.mag();
0114     return v.z() / v.mag();
0115   } else {
0116     LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax << " Signal not contained in strip. Correction is "
0117                                      << thickness / pitch * v.x() / v_xy * v.z() / v.mag() << " instead of "
0118                                      << v.z() / v.mag();
0119     return thickness / pitch * v.x() / v_xy * v.z() / v.mag();
0120   }
0121 }