CMS 3D CMS Logo

SiStripFineDelayTLA Class Reference

#include <DQM/SiStripCommissioningSources/plugins/tracking/SiStripFineDelayTLA.h>

List of all members.

Public Member Functions

std::vector< std::pair
< std::pair< DetId, LocalPoint >,
float > > 
findtrackangle (const Trajectory &traj)
std::vector< std::pair
< std::pair< DetId, LocalPoint >,
float > > 
findtrackangle (const std::vector< Trajectory > &traj)
std::vector< std::pair
< std::pair< DetId, LocalPoint >,
float > > 
findtrackangle (const reco::Track &theT)
std::vector< std::pair
< std::pair< DetId, LocalPoint >,
float > > 
findtrackangle (const TrajectorySeed &seed, const reco::Track &theT)
void init (const edm::Event &e, const edm::EventSetup &c)
 SiStripFineDelayTLA (const edm::ParameterSet &conf)
virtual ~SiStripFineDelayTLA ()

Private Member Functions

double computeAngleCorr (const LocalVector &v, double pitch, double thickness)

Private Attributes

edm::ParameterSet conf_
SimpleTrackRefitterrefitter_
const TrackerGeometrytracker


Detailed Description

Definition at line 21 of file SiStripFineDelayTLA.h.


Constructor & Destructor Documentation

SiStripFineDelayTLA::SiStripFineDelayTLA ( const edm::ParameterSet conf  )  [explicit]

Definition at line 27 of file SiStripFineDelayTLA.cc.

References refitter_.

00027                                                                     : 
00028   conf_(conf)
00029 {
00030   refitter_ = new SimpleTrackRefitter(conf);
00031 }

SiStripFineDelayTLA::~SiStripFineDelayTLA (  )  [virtual]

Definition at line 44 of file SiStripFineDelayTLA.cc.

References refitter_.

00045 {  
00046   delete refitter_;
00047 }  


Member Function Documentation

double SiStripFineDelayTLA::computeAngleCorr ( const LocalVector v,
double  pitch,
double  thickness 
) [private]

Definition at line 128 of file SiStripFineDelayTLA.cc.

References L, LogDebug, PV3DBase< T, PVType, FrameType >::mag(), funct::sqrt(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by findtrackangle().

00129 {
00130   double v_xy = sqrt(v.x()*v.x()+v.y()*v.y());
00131   double L = fabs(thickness*v_xy/v.z());
00132   double Lmax = fabs(pitch/v.x()*v_xy);
00133   if(L<Lmax) {
00134     LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax 
00135         << " Signal contained in strip. Correction is " << v.z()/v.mag();
00136     return v.z()/v.mag();
00137   } else {
00138     LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax
00139        << " Signal not contained in strip. Correction is " << thickness/pitch*v.x()/v_xy*v.z()/v.mag()
00140        << " instead of " << v.z()/v.mag();
00141     return thickness/pitch*v.x()/v_xy*v.z()/v.mag();
00142   }
00143 }

std::vector< std::pair< std::pair< DetId, LocalPoint >,float > > SiStripFineDelayTLA::findtrackangle ( const Trajectory traj  ) 

Definition at line 68 of file SiStripFineDelayTLA.cc.

References angle(), computeAngleCorr(), TrackingRecHit::geographicalId(), TrackerGeometry::idToDet(), TrajectoryStateOnSurface::localDirection(), BaseSiTrackerRecHit2DLocalPos::localPosition(), TrajectoryStateOnSurface::localPosition(), Trajectory::measurements(), GluedGeomDet::monoDet(), SiStripMatchedRecHit2D::monoHit(), GluedGeomDet::stereoDet(), SiStripMatchedRecHit2D::stereoHit(), GeomDet::toGlobal(), GeomDet::toLocal(), GeomDetUnit::topology(), tracker, and PV3DBase< T, PVType, FrameType >::z().

00069 {
00070   std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> >hitangleassociation;
00071   std::vector<TrajectoryMeasurement> TMeas=traj.measurements();
00072   std::vector<TrajectoryMeasurement>::iterator itm;
00073   for (itm=TMeas.begin();itm!=TMeas.end();itm++){
00074     TrajectoryStateOnSurface tsos=itm->updatedState();
00075     const TransientTrackingRecHit::ConstRecHitPointer thit=itm->recHit();
00076     const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>((*thit).hit());
00077     const SiStripRecHit2D* hit=dynamic_cast<const SiStripRecHit2D*>((*thit).hit());
00078     LocalVector trackdirection=tsos.localDirection();
00079     if(matchedhit){//if matched hit...
00080         GluedGeomDet * gdet=(GluedGeomDet *)tracker->idToDet(matchedhit->geographicalId());
00081         GlobalVector gtrkdir=gdet->toGlobal(trackdirection);
00082         // trackdirection on mono det
00083         // this the pointer to the mono hit of a matched hit 
00084         const SiStripRecHit2D *monohit=matchedhit->monoHit();
00085         const GeomDetUnit * monodet=gdet->monoDet();
00086         LocalVector monotkdir=monodet->toLocal(gtrkdir);
00087         if(monotkdir.z()!=0){
00088           // the local angle (mono)
00089           float localpitch = ((StripTopology*)(&monodet->topology()))->localPitch(tsos.localPosition());
00090           float thickness = ((((((monohit->geographicalId())>>25)&0x7f)==0xd)||
00091                              ((((monohit->geographicalId())>>25)&0x7f)==0xe))&&
00092                                    ((((monohit->geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
00093           float angle = computeAngleCorr(monotkdir, localpitch, thickness);
00094           hitangleassociation.push_back(make_pair(make_pair(monohit->geographicalId(),monohit->localPosition()), angle)); 
00095           // trackdirection on stereo det
00096           // this the pointer to the stereo hit of a matched hit 
00097           const SiStripRecHit2D *stereohit=matchedhit->stereoHit();
00098           const GeomDetUnit * stereodet=gdet->stereoDet(); 
00099           LocalVector stereotkdir=stereodet->toLocal(gtrkdir);
00100           if(stereotkdir.z()!=0){
00101             // the local angle (stereo)
00102             float localpitch = ((StripTopology*)(&stereodet->topology()))->localPitch(tsos.localPosition());
00103             float thickness = ((((((stereohit->geographicalId())>>25)&0x7f)==0xd)||
00104                                ((((stereohit->geographicalId())>>25)&0x7f)==0xe))&&
00105                                      ((((stereohit->geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
00106             float angle = computeAngleCorr(stereotkdir, localpitch, thickness);
00107             hitangleassociation.push_back(make_pair(make_pair(stereohit->geographicalId(),stereohit->localPosition()), angle)); 
00108           }
00109         }
00110     }
00111     else if(hit){
00112         GeomDetUnit * det=(GeomDetUnit *)tracker->idToDet(hit->geographicalId());
00113         //  hit= pointer to the rechit
00114         if(trackdirection.z()!=0){
00115           // the local angle (single hit)
00116           float localpitch = ((StripTopology*)(&det->topology()))->localPitch(tsos.localPosition());
00117           float thickness = ((((((hit->geographicalId())>>25)&0x7f)==0xd)||
00118                              ((((hit->geographicalId())>>25)&0x7f)==0xe))&&
00119                                    ((((hit->geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
00120           float angle = computeAngleCorr(trackdirection, localpitch, thickness);
00121           hitangleassociation.push_back(make_pair(make_pair(hit->geographicalId(),hit->localPosition()), angle)); 
00122         }
00123     }
00124   }
00125   return hitangleassociation;
00126 }

std::vector< std::pair< std::pair< DetId, LocalPoint >,float > > SiStripFineDelayTLA::findtrackangle ( const std::vector< Trajectory > &  traj  ) 

Definition at line 60 of file SiStripFineDelayTLA.cc.

References findtrackangle().

00061 {
00062   if (trajVec.size()) {
00063   return findtrackangle(trajVec.front()); }
00064   std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > hitangleassociation;
00065   return hitangleassociation;
00066 }

std::vector< std::pair< std::pair< DetId, LocalPoint >,float > > SiStripFineDelayTLA::findtrackangle ( const reco::Track theT  ) 

Definition at line 49 of file SiStripFineDelayTLA.cc.

References findtrackangle(), refitter_, and SimpleTrackRefitter::refitTrack().

00050 {
00051   std::vector<Trajectory> traj = refitter_->refitTrack(theT);
00052   return findtrackangle(traj);
00053 }

std::vector< std::pair< std::pair< DetId, LocalPoint >,float > > SiStripFineDelayTLA::findtrackangle ( const TrajectorySeed seed,
const reco::Track theT 
)

Definition at line 55 of file SiStripFineDelayTLA.cc.

References refitter_, and SimpleTrackRefitter::refitTrack().

Referenced by SiStripFineDelayHit::detId(), and findtrackangle().

00055                                                                                                                                             {
00056   vector<Trajectory> traj = refitter_->refitTrack(seed,theT);
00057   return findtrackangle(traj);
00058 }       

void SiStripFineDelayTLA::init ( const edm::Event e,
const edm::EventSetup c 
)

Definition at line 33 of file SiStripFineDelayTLA.cc.

References edm::EventSetup::get(), refitter_, SimpleTrackRefitter::setServices(), and tracker.

Referenced by SiStripFineDelayHit::produce().

00034 {
00035   // get geometry
00036   edm::ESHandle<TrackerGeometry> estracker;
00037   es.get<TrackerDigiGeometryRecord>().get(estracker);
00038   tracker=&(*estracker);
00039   // the refitter
00040   refitter_->setServices(es);
00041 }


Member Data Documentation

edm::ParameterSet SiStripFineDelayTLA::conf_ [private]

Definition at line 39 of file SiStripFineDelayTLA.h.

SimpleTrackRefitter* SiStripFineDelayTLA::refitter_ [private]

Definition at line 40 of file SiStripFineDelayTLA.h.

Referenced by findtrackangle(), init(), SiStripFineDelayTLA(), and ~SiStripFineDelayTLA().

const TrackerGeometry* SiStripFineDelayTLA::tracker [private]

Definition at line 41 of file SiStripFineDelayTLA.h.

Referenced by findtrackangle(), and init().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:32:23 2009 for CMSSW by  doxygen 1.5.4