CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

SiStripFineDelayTLA Class Reference

#include <SiStripFineDelayTLA.h>

List of all members.

Public Member Functions

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 Trajectory &traj)
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_
const TrackerGeometrytracker

Detailed Description

Definition at line 20 of file SiStripFineDelayTLA.h.


Constructor & Destructor Documentation

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

Definition at line 26 of file SiStripFineDelayTLA.cc.

                                                                    : 
  conf_(conf)
{
}
SiStripFineDelayTLA::~SiStripFineDelayTLA ( ) [virtual]

Definition at line 40 of file SiStripFineDelayTLA.cc.

{  
}  

Member Function Documentation

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

Definition at line 112 of file SiStripFineDelayTLA.cc.

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

Referenced by findtrackangle().

{
  double v_xy = sqrt(v.x()*v.x()+v.y()*v.y());
  double L = fabs(thickness*v_xy/v.z());
  double Lmax = fabs(pitch/v.x()*v_xy);
  if(L<Lmax) {
    LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax 
        << " Signal contained in strip. Correction is " << v.z()/v.mag();
    return v.z()/v.mag();
  } else {
    LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax
       << " Signal not contained in strip. Correction is " << thickness/pitch*v.x()/v_xy*v.z()/v.mag()
       << " instead of " << v.z()/v.mag();
    return thickness/pitch*v.x()/v_xy*v.z()/v.mag();
  }
}
std::vector< std::pair< std::pair< DetId, LocalPoint >,float > > SiStripFineDelayTLA::findtrackangle ( const std::vector< Trajectory > &  traj)

Definition at line 44 of file SiStripFineDelayTLA.cc.

Referenced by SiStripFineDelayHit::detId().

{
  if (trajVec.size()) {
  return findtrackangle(trajVec.front()); }
  std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > hitangleassociation;
  return hitangleassociation;
}
std::vector< std::pair< std::pair< DetId, LocalPoint >,float > > SiStripFineDelayTLA::findtrackangle ( const Trajectory traj)

Definition at line 52 of file SiStripFineDelayTLA.cc.

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

{
  std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> >hitangleassociation;
  std::vector<TrajectoryMeasurement> TMeas=traj.measurements();
  std::vector<TrajectoryMeasurement>::iterator itm;
  for (itm=TMeas.begin();itm!=TMeas.end();itm++){
    TrajectoryStateOnSurface tsos=itm->updatedState();
    const TransientTrackingRecHit::ConstRecHitPointer thit=itm->recHit();
    const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>((*thit).hit());
    const SiStripRecHit2D* hit=dynamic_cast<const SiStripRecHit2D*>((*thit).hit());
    LocalVector trackdirection=tsos.localDirection();
    if(matchedhit){//if matched hit...
        GluedGeomDet * gdet=(GluedGeomDet *)tracker->idToDet(matchedhit->geographicalId());
        GlobalVector gtrkdir=gdet->toGlobal(trackdirection);
        // trackdirection on mono det
        // this the pointer to the mono hit of a matched hit 
        const SiStripRecHit2D monohit=matchedhit->monoHit();
        const GeomDetUnit * monodet=gdet->monoDet();
        LocalVector monotkdir=monodet->toLocal(gtrkdir);
        if(monotkdir.z()!=0){
          // the local angle (mono)
          float localpitch = ((StripTopology*)(&monodet->topology()))->localPitch(tsos.localPosition());
          float thickness = ((((((monohit.geographicalId())>>25)&0x7f)==0xd)||
                             ((((monohit.geographicalId())>>25)&0x7f)==0xe))&&
                                   ((((monohit.geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
          float angle = computeAngleCorr(monotkdir, localpitch, thickness);
          hitangleassociation.push_back(make_pair(make_pair(monohit.geographicalId(),monohit.localPosition()), angle)); 
          // trackdirection on stereo det
          // this the pointer to the stereo hit of a matched hit 
          const SiStripRecHit2D stereohit=matchedhit->stereoHit();
          const GeomDetUnit * stereodet=gdet->stereoDet(); 
          LocalVector stereotkdir=stereodet->toLocal(gtrkdir);
          if(stereotkdir.z()!=0){
            // the local angle (stereo)
            float localpitch = ((StripTopology*)(&stereodet->topology()))->localPitch(tsos.localPosition());
            float thickness = ((((((stereohit.geographicalId())>>25)&0x7f)==0xd)||
                               ((((stereohit.geographicalId())>>25)&0x7f)==0xe))&&
                                     ((((stereohit.geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
            float angle = computeAngleCorr(stereotkdir, localpitch, thickness);
            hitangleassociation.push_back(make_pair(make_pair(stereohit.geographicalId(),stereohit.localPosition()), angle)); 
          }
        }
    }
    else if(hit){
        GeomDetUnit * det=(GeomDetUnit *)tracker->idToDet(hit->geographicalId());
        //  hit= pointer to the rechit
        if(trackdirection.z()!=0){
          // the local angle (single hit)
          float localpitch = ((StripTopology*)(&det->topology()))->localPitch(tsos.localPosition());
          float thickness = ((((((hit->geographicalId())>>25)&0x7f)==0xd)||
                             ((((hit->geographicalId())>>25)&0x7f)==0xe))&&
                                   ((((hit->geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
          float angle = computeAngleCorr(trackdirection, localpitch, thickness);
          hitangleassociation.push_back(make_pair(make_pair(hit->geographicalId(),hit->localPosition()), angle)); 
        }
    }
  }
  return hitangleassociation;
}
void SiStripFineDelayTLA::init ( const edm::Event e,
const edm::EventSetup c 
)

Definition at line 31 of file SiStripFineDelayTLA.cc.

References edm::EventSetup::get(), and tracker.

Referenced by SiStripFineDelayHit::produce().

{
  // get geometry
  edm::ESHandle<TrackerGeometry> estracker;
  es.get<TrackerDigiGeometryRecord>().get(estracker);
  tracker=&(*estracker);
}

Member Data Documentation

Definition at line 36 of file SiStripFineDelayTLA.h.

Definition at line 37 of file SiStripFineDelayTLA.h.

Referenced by findtrackangle(), and init().