#include <IdealHelixParameters.h>
Public Member Functions | |
math::XYZVector | GetCircleCenter () const |
math::XYZVector | GetMomentumAtTangentPoint () const |
float | GetRotationAngle () const |
math::XYZVector | GetTangentPoint () const |
float | GetTransverseIPAtTangent () const |
IdealHelixParameters () | |
bool | isTangentPointDistanceLessThan (float rmax, const reco::Track *track, const math::XYZVector refPoint) |
void | setData (const reco::Track *track, math::XYZVector refPoint=math::XYZVector(0, 0, 0)) |
void | setData (const reco::Track *track, math::XYZPoint ref) |
void | setMagnField (const MagneticField *magnField) |
~IdealHelixParameters () | |
Private Member Functions | |
void | calculate () |
void | evalCircleCenter () |
void | evalMomentumatTangentPoint () |
void | evalTangentPoint () |
Private Attributes | |
math::XYZVector | _circleCenter |
const MagneticField * | _magnField |
math::XYZVector | _MomentumAtTangentPoint |
float | _radius |
math::XYZVector | _refPoint |
float | _rotationAngle |
math::XYZVector | _tangentPoint |
const reco::Track * | _track |
float | _transverseIP |
Definition at line 21 of file IdealHelixParameters.h.
IdealHelixParameters::IdealHelixParameters | ( | ) | [inline] |
Definition at line 25 of file IdealHelixParameters.h.
: _magnField(0),_track(0), _refPoint (math::XYZVector(0,0,0)), _tangentPoint (math::XYZVector(0,0,0)), _MomentumAtTangentPoint(math::XYZVector(0,0,0)){};
IdealHelixParameters::~IdealHelixParameters | ( | ) | [inline] |
Definition at line 30 of file IdealHelixParameters.h.
{};
void IdealHelixParameters::calculate | ( | ) | [private] |
Definition at line 25 of file IdealHelixParameters.cc.
References _refPoint, _track, evalCircleCenter(), evalMomentumatTangentPoint(), evalTangentPoint(), reco::Track::innerDetId(), reco::Track::innerMomentum(), reco::Track::innerPosition(), and LogDebug.
Referenced by setData().
{ LogDebug("IdealHelixParameters") << "[IdealHelixParameters] " << "\n\t "<< "refPoint \t" <<_refPoint << "\n\t "<< "innerDetid \t" << _track->innerDetId() << "\n\t "<< "innerMomentum \t" << _track->innerMomentum() << "\n\t "<< "innerPosition \t" << _track->innerPosition(); evalCircleCenter(); evalTangentPoint(); evalMomentumatTangentPoint(); /* DEBUG math::XYZVector _checkMomentumAtTangentoPoint ( cos_angle*innerMomentum.x()+sin_angle*innerMomentum.y(), -1*sin_angle*innerMomentum.x()+cos_angle*innerMomentum.y(), innerMomentum.z() ); math::XYZVector cp_xy(_checkMomentumAtTangentoPoint.x(),_checkMomentumAtTangentoPoint.y(),0); float _checktransverseIP = (r_xy.Cross(cp_xy)).z()/cp_xy.rho(); ss << "\n\t "<< "_checkMomentumAtTangentoPoint \t" <<_checkMomentumAtTangentoPoint << "\n\t "<< "_checktransverseIP \t" <<_checktransverseIP; */ }
void IdealHelixParameters::evalCircleCenter | ( | ) | [private] |
Definition at line 58 of file IdealHelixParameters.cc.
References _circleCenter, _magnField, _radius, _track, reco::TrackBase::charge(), funct::cos(), reco::Track::innerMomentum(), reco::Track::innerPosition(), LogDebug, phi, and funct::sin().
Referenced by calculate().
{ GlobalPoint innerPos(_track->innerPosition().x(),_track->innerPosition().y(),_track->innerPosition().z()); GlobalVector innerMom(_track->innerMomentum().x(),_track->innerMomentum().y(),_track->innerMomentum().z()); GlobalTrajectoryParameters globalTrajParam( innerPos, innerMom, _track->charge(), _magnField ); _radius = 1./fabs(globalTrajParam.transverseCurvature()); float phi = innerMom.phi(); //eval Center of the Circumference passing in track->innserPosition float Xc = innerPos.x() + _track->charge()* _radius * sin(phi); float Yc = innerPos.y() - _track->charge()* _radius * cos(phi); _circleCenter.SetXYZ(Xc,Yc,innerPos.z()); //NB Z component is useless LogDebug("IdealHelixParameters") << "\n\t "<< "circle center \t" <<_circleCenter << "\n\t "<< "radius " <<_radius; }
void IdealHelixParameters::evalMomentumatTangentPoint | ( | ) | [private] |
Definition at line 87 of file IdealHelixParameters.cc.
References _circleCenter, _MomentumAtTangentPoint, _rotationAngle, _tangentPoint, _track, _transverseIP, reco::Track::innerMomentum(), reco::Track::innerPosition(), and LogDebug.
Referenced by calculate().
{ if(_tangentPoint.r()==0){ _MomentumAtTangentPoint = math::XYZVector(0.,0.,0.); return; } math::XYZVector innerPosition(_track->innerPosition().X(),_track->innerPosition().Y(),_track->innerPosition().Z()); math::XYZVector innerMomentum(_track->innerMomentum()); math::XYZVector pCi=innerPosition-_circleCenter; math::XYZVector pCT=_tangentPoint-_circleCenter; math::XYZVector pCi_xy(pCi.x(),pCi.y(),0); math::XYZVector pCT_xy(pCT.x(),pCT.y(),0); float cos_angle = pCi_xy.Dot(pCT_xy)/pCi_xy.rho()/pCT_xy.rho(); float sin_angle = pCi_xy.Cross(pCT_xy).z()/pCi_xy.rho()/pCT_xy.rho(); _MomentumAtTangentPoint = math::XYZVector( cos_angle*innerMomentum.x()-sin_angle*innerMomentum.y(), sin_angle*innerMomentum.x()+cos_angle*innerMomentum.y(), innerMomentum.z() ); math::XYZVector r_(_tangentPoint.x(),_tangentPoint.y(),0); math::XYZVector p_(_MomentumAtTangentPoint.x(),_MomentumAtTangentPoint.y(),0); _transverseIP = (r_.x()*p_.y()-r_.y()*p_.x())/p_.rho()/r_.rho(); _rotationAngle=atan(sin_angle/cos_angle); LogDebug("IdealHelixParameters") << "\n\t "<< "_MomentumAtTangentPoint \t" <<_MomentumAtTangentPoint << "\n\t "<< "sin_angle \t" <<sin_angle << " angle \t" <<asin(sin_angle) << "\n\t "<< "cos_angle \t" <<cos_angle << " angle \t" <<acos(cos_angle) << "\n\t "<< "_rotationAngle \t" <<_rotationAngle << "\n\t "<< "_transverseIP \t" <<_transverseIP << "\n\t "<< " check similitude pz/pt "<< _track->innerMomentum().z()/_track->innerMomentum().rho() << " pZ/pRho " << innerPosition.z()/innerPosition.rho(); }
void IdealHelixParameters::evalTangentPoint | ( | ) | [private] |
Definition at line 132 of file IdealHelixParameters.cc.
References _circleCenter, _radius, _refPoint, _tangentPoint, _track, reco::Track::innerMomentum(), reco::Track::innerPosition(), prof2calltree::l, LogDebug, rho, and mathSSE::sqrt().
Referenced by calculate().
{ math::XYZVector innerPosition(_track->innerPosition().X(),_track->innerPosition().Y(),_track->innerPosition().Z()); math::XYZVector vL = _circleCenter-_refPoint; float l = vL.rho(); if(l<_radius){ //refpoint is inside the circle _tangentPoint.SetXYZ(0.,0.,0.); return; } float sin_alpha = _radius/l; //there are two possible tangents, with Point of tangence T1 and T2, and same distance T from _refPoint float T=sqrt(l*l-_radius*_radius); math::XYZVector vLperp(-vL.y(),vL.x(),0); //NB: z component not correct math::XYZVector tmpT1 = (1-sin_alpha*sin_alpha)*vL + T/l*sin_alpha*vLperp; math::XYZVector tmpT2 = (1-sin_alpha*sin_alpha)*vL - T/l*sin_alpha*vLperp; if( (tmpT1-innerPosition).rho()<(tmpT2-innerPosition).rho()) _tangentPoint=tmpT1+_refPoint; else _tangentPoint=tmpT2+_refPoint; //Fix the Z component _tangentPoint.SetZ( (_tangentPoint.rho()-_refPoint.rho()) * (_track->innerMomentum().z()/_track->innerMomentum().rho()) + _refPoint.z() ); LogDebug("IdealHelixParameters")<< "\n\t "<< "tangent Point \t" <<_tangentPoint; }
math::XYZVector IdealHelixParameters::GetCircleCenter | ( | ) | const [inline] |
Definition at line 38 of file IdealHelixParameters.h.
References _circleCenter.
{return _circleCenter;}
math::XYZVector IdealHelixParameters::GetMomentumAtTangentPoint | ( | ) | const [inline] |
Definition at line 40 of file IdealHelixParameters.h.
References _MomentumAtTangentPoint.
Referenced by PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::inspectTrack().
{return _MomentumAtTangentPoint;}
float IdealHelixParameters::GetRotationAngle | ( | ) | const [inline] |
Definition at line 42 of file IdealHelixParameters.h.
References _rotationAngle.
{return _rotationAngle;}
math::XYZVector IdealHelixParameters::GetTangentPoint | ( | ) | const [inline] |
Definition at line 39 of file IdealHelixParameters.h.
References _tangentPoint.
Referenced by PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::inspectTrack(), isTangentPointDistanceLessThan(), and PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::rejectTrack().
{return _tangentPoint;}
float IdealHelixParameters::GetTransverseIPAtTangent | ( | ) | const [inline] |
Definition at line 41 of file IdealHelixParameters.h.
References _transverseIP.
{return _transverseIP;}
bool IdealHelixParameters::isTangentPointDistanceLessThan | ( | float | rmax, |
const reco::Track * | track, | ||
const math::XYZVector | refPoint | ||
) |
Definition at line 166 of file IdealHelixParameters.cc.
References GetTangentPoint(), alignCSCRings::r, rho, and setData().
Referenced by ConversionTrackProducer::produce().
{ setData(track,refPoint); if(GetTangentPoint().r()==0){ //this case means a null results on the IdealHelixParameters side return true; } if(GetTangentPoint().rho()<rMax){ //this case means a track that has the tangent point nearby the primary vertex // if the track is primary, this number tends to be the primary vertex itself //Rejecting all the potential photon conversions having a "vertex" inside the beampipe //We should not miss too much, seen that the conversions at the beam pipe are the better reconstructed return true; } return false; }
void IdealHelixParameters::setData | ( | const reco::Track * | track, |
math::XYZVector | refPoint = math::XYZVector(0,0,0) |
||
) |
Definition at line 9 of file IdealHelixParameters.cc.
References _refPoint, _track, calculate(), and alignCSCRings::r.
Referenced by PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::inspectTrack(), isTangentPointDistanceLessThan(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::rejectTrack(), and setData().
void IdealHelixParameters::setData | ( | const reco::Track * | track, |
math::XYZPoint | ref | ||
) |
Definition at line 4 of file IdealHelixParameters.cc.
References setData().
{ setData(track,math::XYZVector(ref.x(),ref.y(),ref.z())); }
void IdealHelixParameters::setMagnField | ( | const MagneticField * | magnField | ) | [inline] |
Definition at line 32 of file IdealHelixParameters.h.
References _magnField.
Referenced by PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::analyze(), and ConversionTrackProducer::produce().
{_magnField=magnField;}
Definition at line 54 of file IdealHelixParameters.h.
Referenced by evalCircleCenter(), evalMomentumatTangentPoint(), evalTangentPoint(), and GetCircleCenter().
const MagneticField* IdealHelixParameters::_magnField [private] |
Definition at line 51 of file IdealHelixParameters.h.
Referenced by evalCircleCenter(), and setMagnField().
Definition at line 57 of file IdealHelixParameters.h.
Referenced by evalMomentumatTangentPoint(), and GetMomentumAtTangentPoint().
float IdealHelixParameters::_radius [private] |
Definition at line 53 of file IdealHelixParameters.h.
Referenced by evalCircleCenter(), and evalTangentPoint().
Definition at line 55 of file IdealHelixParameters.h.
Referenced by calculate(), evalTangentPoint(), and setData().
float IdealHelixParameters::_rotationAngle [private] |
Definition at line 59 of file IdealHelixParameters.h.
Referenced by evalMomentumatTangentPoint(), and GetRotationAngle().
Definition at line 56 of file IdealHelixParameters.h.
Referenced by evalMomentumatTangentPoint(), evalTangentPoint(), and GetTangentPoint().
const reco::Track* IdealHelixParameters::_track [private] |
Definition at line 52 of file IdealHelixParameters.h.
Referenced by calculate(), evalCircleCenter(), evalMomentumatTangentPoint(), evalTangentPoint(), and setData().
float IdealHelixParameters::_transverseIP [private] |
Definition at line 58 of file IdealHelixParameters.h.
Referenced by evalMomentumatTangentPoint(), and GetTransverseIPAtTangent().