CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTLA.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include <string>
00003 #include <iostream>
00004 #include <TMath.h>
00005 #include "DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTLA.h"
00006 
00007 #include "FWCore/Framework/interface/EventSetup.h"
00008 #include "FWCore/Framework/interface/ESHandle.h"
00009 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00010 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00011 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00012 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00013 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00014 #include "Geometry/TrackerGeometryBuilder/interface/GluedGeomDet.h"
00015 #include <Geometry/CommonTopologies/interface/Topology.h>
00016 #include <Geometry/CommonTopologies/interface/StripTopology.h>
00017 #include "DataFormats/TrackReco/interface/Track.h"
00018 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00019 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00020 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
00021 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00022 #include <TrackingTools/PatternTools/interface/Trajectory.h>
00023 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00024 
00025 using namespace std;
00026 SiStripFineDelayTLA::SiStripFineDelayTLA(edm::ParameterSet const& conf) : 
00027   conf_(conf)
00028 {
00029 }
00030 
00031 void SiStripFineDelayTLA::init(const edm::Event& e, const edm::EventSetup& es)
00032 {
00033   // get geometry
00034   edm::ESHandle<TrackerGeometry> estracker;
00035   es.get<TrackerDigiGeometryRecord>().get(estracker);
00036   tracker=&(*estracker);
00037 }
00038 
00039 // Virtual destructor needed.
00040 SiStripFineDelayTLA::~SiStripFineDelayTLA() 
00041 {  
00042 }  
00043 
00044 std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > SiStripFineDelayTLA::findtrackangle(const std::vector<Trajectory>& trajVec)
00045 {
00046   if (trajVec.size()) {
00047   return findtrackangle(trajVec.front()); }
00048   std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > hitangleassociation;
00049   return hitangleassociation;
00050 }
00051 
00052 std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > SiStripFineDelayTLA::findtrackangle(const Trajectory& traj)
00053 {
00054   std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> >hitangleassociation;
00055   std::vector<TrajectoryMeasurement> TMeas=traj.measurements();
00056   std::vector<TrajectoryMeasurement>::iterator itm;
00057   for (itm=TMeas.begin();itm!=TMeas.end();itm++){
00058     TrajectoryStateOnSurface tsos=itm->updatedState();
00059     const TransientTrackingRecHit::ConstRecHitPointer thit=itm->recHit();
00060     const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>((*thit).hit());
00061     const SiStripRecHit2D* hit=dynamic_cast<const SiStripRecHit2D*>((*thit).hit());
00062     LocalVector trackdirection=tsos.localDirection();
00063     if(matchedhit){//if matched hit...
00064         GluedGeomDet * gdet=(GluedGeomDet *)tracker->idToDet(matchedhit->geographicalId());
00065         GlobalVector gtrkdir=gdet->toGlobal(trackdirection);
00066         // trackdirection on mono det
00067         // this the pointer to the mono hit of a matched hit 
00068         const SiStripRecHit2D monohit=matchedhit->monoHit();
00069         const GeomDetUnit * monodet=gdet->monoDet();
00070         LocalVector monotkdir=monodet->toLocal(gtrkdir);
00071         if(monotkdir.z()!=0){
00072           // the local angle (mono)
00073           float localpitch = ((StripTopology*)(&monodet->topology()))->localPitch(tsos.localPosition());
00074           float thickness = ((((((monohit.geographicalId())>>25)&0x7f)==0xd)||
00075                              ((((monohit.geographicalId())>>25)&0x7f)==0xe))&&
00076                                    ((((monohit.geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
00077           float angle = computeAngleCorr(monotkdir, localpitch, thickness);
00078           hitangleassociation.push_back(make_pair(make_pair(monohit.geographicalId(),monohit.localPosition()), angle)); 
00079           // trackdirection on stereo det
00080           // this the pointer to the stereo hit of a matched hit 
00081           const SiStripRecHit2D stereohit=matchedhit->stereoHit();
00082           const GeomDetUnit * stereodet=gdet->stereoDet(); 
00083           LocalVector stereotkdir=stereodet->toLocal(gtrkdir);
00084           if(stereotkdir.z()!=0){
00085             // the local angle (stereo)
00086             float localpitch = ((StripTopology*)(&stereodet->topology()))->localPitch(tsos.localPosition());
00087             float thickness = ((((((stereohit.geographicalId())>>25)&0x7f)==0xd)||
00088                                ((((stereohit.geographicalId())>>25)&0x7f)==0xe))&&
00089                                      ((((stereohit.geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
00090             float angle = computeAngleCorr(stereotkdir, localpitch, thickness);
00091             hitangleassociation.push_back(make_pair(make_pair(stereohit.geographicalId(),stereohit.localPosition()), angle)); 
00092           }
00093         }
00094     }
00095     else if(hit){
00096         GeomDetUnit * det=(GeomDetUnit *)tracker->idToDet(hit->geographicalId());
00097         //  hit= pointer to the rechit
00098         if(trackdirection.z()!=0){
00099           // the local angle (single hit)
00100           float localpitch = ((StripTopology*)(&det->topology()))->localPitch(tsos.localPosition());
00101           float thickness = ((((((hit->geographicalId())>>25)&0x7f)==0xd)||
00102                              ((((hit->geographicalId())>>25)&0x7f)==0xe))&&
00103                                    ((((hit->geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
00104           float angle = computeAngleCorr(trackdirection, localpitch, thickness);
00105           hitangleassociation.push_back(make_pair(make_pair(hit->geographicalId(),hit->localPosition()), angle)); 
00106         }
00107     }
00108   }
00109   return hitangleassociation;
00110 }
00111 
00112 double SiStripFineDelayTLA::computeAngleCorr(const LocalVector& v, double pitch, double thickness)
00113 {
00114   double v_xy = sqrt(v.x()*v.x()+v.y()*v.y());
00115   double L = fabs(thickness*v_xy/v.z());
00116   double Lmax = fabs(pitch/v.x()*v_xy);
00117   if(L<Lmax) {
00118     LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax 
00119         << " Signal contained in strip. Correction is " << v.z()/v.mag();
00120     return v.z()/v.mag();
00121   } else {
00122     LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax
00123        << " Signal not contained in strip. Correction is " << thickness/pitch*v.x()/v_xy*v.z()/v.mag()
00124        << " instead of " << v.z()/v.mag();
00125     return thickness/pitch*v.x()/v_xy*v.z()/v.mag();
00126   }
00127 }
00128