CMS 3D CMS Logo

StripCPEfromTrackAngle Class Reference

#include <RecoLocalTracker/SiStripRecHitConverter/interface/StripCPEfromTrackAngle.h>

Inheritance diagram for StripCPEfromTrackAngle:

StripCPE StripClusterParameterEstimator ClusterParameterEstimator< SiStripCluster >

List of all members.

Public Member Functions

StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &cl, const LocalTrajectoryParameters &ltp) const
StripClusterParameterEstimator::LocalValues localParameters (const SiStripCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
 StripCPEfromTrackAngle (edm::ParameterSet &conf, const MagneticField *mag, const TrackerGeometry *geom, const SiStripLorentzAngle *LorentzAngle)


Detailed Description

Definition at line 12 of file StripCPEfromTrackAngle.h.


Constructor & Destructor Documentation

StripCPEfromTrackAngle::StripCPEfromTrackAngle ( edm::ParameterSet conf,
const MagneticField mag,
const TrackerGeometry geom,
const SiStripLorentzAngle LorentzAngle 
) [inline]

Definition at line 16 of file StripCPEfromTrackAngle.h.

00020     :StripCPE(conf,mag, geom, LorentzAngle ){}


Member Function Documentation

StripClusterParameterEstimator::LocalValues StripCPEfromTrackAngle::localParameters ( const SiStripCluster cl,
const LocalTrajectoryParameters ltp 
) const

Definition at line 6 of file StripCPEfromTrackAngle.cc.

References SiStripCluster::amplitudes(), SiStripCluster::barycenter(), StripCPE::Param::drift, funct::exp(), SiStripCluster::geographicalId(), StripTopology::localError(), StripTopology::localPosition(), m1, PV3DBase< T, PVType, FrameType >::mag(), StripCPE::Param::maxLength, Topology::measurementPosition(), min, LocalTrajectoryParameters::momentum(), StripCPE::Param::nstrips, p, p1, p2, StripCPE::param(), LocalTrajectoryParameters::position(), HLT_VtxMuL3::result, StripCPE::Param::thickness, StripCPE::Param::topology, PV3DBase< T, PVType, FrameType >::x(), PV2DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

00006                                                                                                                                                       {
00007   //
00008   // get the det from the geometry
00009   //
00010 
00011   StripCPE::Param const & p = param(DetId(cl.geographicalId()));
00012   
00013 
00014  
00015   // ionisation length
00016    
00017   //float ionLen = std::min( trackDir.mag(), maxLength);
00018 
00019   //  const float par1=38.07;
00020   // const float par2=0.3184; 
00021   // const float par3=0.09828; 
00022   // const float P1 = par1 * p.thickness; 
00023   // const float P2 = par2; 
00024   // const float P3 = par3;
00025 
00026 
00027   //  drift*=(thickness/2);
00028 
00029   //calculate error form track angle
00030 
00031 
00032   LocalPoint  middlepoint = ltp.position();
00033   LocalVector trackDir = ltp.momentum()/ltp.momentum().mag();
00034   LocalPoint  position = p.topology->localPosition(cl.barycenter());
00035 
00036 
00037 
00038   if(trackDir.z()*p.drift.z() > 0.) trackDir *= -1.;
00039 
00040 
00041   if(trackDir.z() !=0.) {
00042     trackDir *= fabs(p.thickness/trackDir.z());
00043   } else {
00044     trackDir *= p.maxLength/trackDir.mag();
00045   }
00046 
00047 
00048       
00049   if(p.drift.z() == 0.) {
00050     //  if(drift.z() == 0.||cl.amplitudes().size()==1) {
00051     edm::LogError("StripCPE") <<"No drift towards anodes !!!";
00052     LocalError eresult = p.topology->localError(cl.barycenter(),1/12.);
00053     //  LocalPoint  result=LocalPoint(position.x()-drift.x()/2,position.y()-drift.y()/2,0);
00054     return std::make_pair(position-p.drift*(p.thickness/2),eresult);
00055   }      
00056 
00057 
00058 
00059 
00060   // covered length along U
00061       
00062   LocalVector middleOfProjection = 0.5*(trackDir + p.drift);
00063 
00064   LocalPoint middlePointOnStrips = middlepoint + 0.5*p.drift;
00065 
00066   LocalPoint p1 = LocalPoint(middlePointOnStrips.x() + middleOfProjection.x()
00067                              ,middlePointOnStrips.y() + middleOfProjection.y());
00068   LocalPoint p2 = LocalPoint(middlePointOnStrips.x() - middleOfProjection.x()
00069                              ,middlePointOnStrips.y() - middleOfProjection.y());
00070 
00071   MeasurementPoint m1 = p.topology->measurementPosition(p1);
00072   MeasurementPoint m2 = p.topology->measurementPosition(p2);
00073   float u1 = m1.x();
00074   float u2 = m2.x();
00075   float uerr2;
00076   float uProj = std::min( float(fabs( u1 - u2)), float(p.nstrips));
00077   if((cl.amplitudes().size() - (uProj+2.5)) > 1)uerr2=cl.amplitudes().size() * cl.amplitudes().size() * (1./12.);
00078   else{
00079     //   float par1=38.07;
00080     //   float par2=0.3184; 
00081     //   float par3=0.09828; 
00082     //   float P1 = par1 * thickness; 
00083     //   float P2 = par2; 
00084     //   float P3 = par3;
00085     //   float uerr;
00086     
00087     //   uerr =(uProj-P1)*(uProj-P1)*(P2-P3)/(P1*P1)+P3;
00088     
00089     const float P1=-0.339;
00090     const float P2=0.90;
00091     const float P3=0.279;
00092     
00093     float uerr=P1*uProj*exp(-uProj*P2)+P3;
00094     uerr2=uerr*uerr;
00095   }
00096 
00097   MeasurementError merror=MeasurementError( uerr2, 0., 1./12.);
00098   LocalPoint result=LocalPoint(position.x()-0.5*p.drift.x(),position.y()-0.5*p.drift.y(),0);
00099   MeasurementPoint mpoint=p.topology->measurementPosition(result);
00100   LocalError eresult=p.topology->localError(mpoint,merror);
00101   return std::make_pair(result,eresult);
00102 }

StripClusterParameterEstimator::LocalValues StripCPEfromTrackAngle::localParameters ( const SiStripCluster cl,
const GeomDetUnit det,
const LocalTrajectoryParameters ltp 
) const [inline, virtual]

Reimplemented from ClusterParameterEstimator< SiStripCluster >.

Definition at line 25 of file StripCPEfromTrackAngle.h.

00027                                                                                                            {
00028     return localParameters(cl,ltp);
00029   }; 


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