CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
TangentApproachInRPhi Class Reference

#include <TangentApproachInRPhi.h>

Inheritance diagram for TangentApproachInRPhi:
ClosestApproachOnHelices

Public Member Functions

virtual bool calculate (const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
 
virtual bool calculate (const FreeTrajectoryState &sta, const FreeTrajectoryState &stb)
 
virtual TangentApproachInRPhiclone () const
 
virtual GlobalPoint crossingPoint () const
 
virtual float distance () const
 
float perpdist () const
 
virtual std::pair< GlobalPoint,
GlobalPoint
points () const
 
virtual bool status () const
 
 TangentApproachInRPhi ()
 
std::pair
< GlobalTrajectoryParameters,
GlobalTrajectoryParameters
trajectoryParameters () const
 
- Public Member Functions inherited from ClosestApproachOnHelices
 ClosestApproachOnHelices ()
 
virtual ~ClosestApproachOnHelices ()
 

Private Member Functions

bool calculate (const TrackCharge &chargeA, const GlobalVector &momentumA, const GlobalPoint &positionA, const TrackCharge &chargeB, const GlobalVector &momentumB, const GlobalPoint &positionB, const MagneticField &magField)
 
void circleParameters (const TrackCharge &charge, const GlobalVector &momemtum, const GlobalPoint &position, double &xc, double &yc, double &r, const MagneticField &magField) const
 
GlobalTrajectoryParameters trajectoryParameters (const GlobalPoint &newpt, const GlobalTrajectoryParameters &oldpar) const
 
int transverseCoord (double cxa, double cya, double ra, double cxb, double cyb, double rb, double &xg1, double &yg1, double &xg2, double &yg2) const
 
double zCoord (const GlobalVector &mom, const GlobalPoint &pos, double r, double xc, double yc, double xg, double yg) const
 

Private Attributes

bool intersection_
 
GlobalTrajectoryParameters paramA
 
GlobalTrajectoryParameters paramB
 
GlobalPoint posA
 
GlobalPoint posB
 
bool status_
 

Detailed Description

Given two trajectory states, computes the two points of tangent approach in the transverse plane for the helices extrapolated from these states. 1) computes the tangent points of the circles in transverse plane. Two cases: - Whether or not the circles intersect the points used are the points of tangent approach of the two circles. 2) computes the corresponding z-coordinates.

Definition at line 15 of file TangentApproachInRPhi.h.

Constructor & Destructor Documentation

TangentApproachInRPhi::TangentApproachInRPhi ( )
inline

Definition at line 19 of file TangentApproachInRPhi.h.

References intersection_, and status_.

Referenced by clone().

Member Function Documentation

bool TangentApproachInRPhi::calculate ( const TrajectoryStateOnSurface sta,
const TrajectoryStateOnSurface stb 
)
virtual

Implements ClosestApproachOnHelices.

Definition at line 7 of file TangentApproachInRPhi.cc.

References TrajectoryStateOnSurface::charge(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalParameters(), TrajectoryStateOnSurface::globalPosition(), GlobalTrajectoryParameters::magneticField(), and FreeTrajectoryState::parameters().

Referenced by ConversionProducer::preselectTrackPair().

9 {
10  TrackCharge chargeA = sta.charge(); TrackCharge chargeB = stb.charge();
11  GlobalVector momentumA = sta.globalMomentum();
12  GlobalVector momentumB = stb.globalMomentum();
13  GlobalPoint positionA = sta.globalPosition();
14  GlobalPoint positionB = stb.globalPosition();
15  paramA = sta.globalParameters();
16  paramB = stb.globalParameters();
17 
18  return calculate(chargeA, momentumA, positionA, chargeB, momentumB, positionB,
20 }
const GlobalTrajectoryParameters & parameters() const
GlobalTrajectoryParameters paramA
GlobalPoint globalPosition() const
int TrackCharge
Definition: TrackCharge.h:4
GlobalTrajectoryParameters paramB
FreeTrajectoryState const * freeState(bool withErrors=true) const
const GlobalTrajectoryParameters & globalParameters() const
GlobalVector globalMomentum() const
const MagneticField & magneticField() const
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
bool TangentApproachInRPhi::calculate ( const FreeTrajectoryState sta,
const FreeTrajectoryState stb 
)
virtual

Implements ClosestApproachOnHelices.

Definition at line 23 of file TangentApproachInRPhi.cc.

References FreeTrajectoryState::charge(), GlobalTrajectoryParameters::magneticField(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), and FreeTrajectoryState::position().

25 {
26  TrackCharge chargeA = sta.charge(); TrackCharge chargeB = stb.charge();
27  GlobalVector momentumA = sta.momentum();
28  GlobalVector momentumB = stb.momentum();
29  GlobalPoint positionA = sta.position();
30  GlobalPoint positionB = stb.position();
31  paramA = sta.parameters();
32  paramB = stb.parameters();
33 
34  return calculate(chargeA, momentumA, positionA, chargeB, momentumB, positionB,
35  sta.parameters().magneticField());
36 }
const GlobalTrajectoryParameters & parameters() const
GlobalTrajectoryParameters paramA
TrackCharge charge() const
int TrackCharge
Definition: TrackCharge.h:4
GlobalTrajectoryParameters paramB
GlobalVector momentum() const
GlobalPoint position() const
const MagneticField & magneticField() const
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
bool TangentApproachInRPhi::calculate ( const TrackCharge chargeA,
const GlobalVector momentumA,
const GlobalPoint positionA,
const TrackCharge chargeB,
const GlobalVector momentumB,
const GlobalPoint positionB,
const MagneticField magField 
)
private

Definition at line 80 of file TangentApproachInRPhi.cc.

87 {
88 
89  // centres and radii of track circles
90  double xca, yca, ra;
91  circleParameters(chargeA, momentumA, positionA, xca, yca, ra, magField);
92  double xcb, ycb, rb;
93  circleParameters(chargeB, momentumB, positionB, xcb, ycb, rb, magField);
94 
95  // points of closest approach in transverse plane
96  double xg1, yg1, xg2, yg2;
97  int flag = transverseCoord(xca, yca, ra, xcb, ycb, rb, xg1, yg1, xg2, yg2);
98  if (flag == 0) {
99  status_ = false;
100  return false;
101  }
102 
103  double xga, yga, zga, xgb, ygb, zgb;
104 
105  if (flag == 1) {
106  intersection_ = true;
107  }
108  else {
109  intersection_ = false;
110  }
111 
112  // one point of closest approach on each track in transverse plane
113  xga = xg1; yga = yg1;
114  zga = zCoord(momentumA, positionA, ra, xca, yca, xga, yga);
115  xgb = xg2; ygb = yg2;
116  zgb = zCoord(momentumB, positionB, rb, xcb, ycb, xgb, ygb);
117 
118 
119  posA = GlobalPoint(xga, yga, zga);
120  posB = GlobalPoint(xgb, ygb, zgb);
121  status_ = true;
122  return true;
123 }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
double zCoord(const GlobalVector &mom, const GlobalPoint &pos, double r, double xc, double yc, double xg, double yg) const
int transverseCoord(double cxa, double cya, double ra, double cxb, double cyb, double rb, double &xg1, double &yg1, double &xg2, double &yg2) const
void circleParameters(const TrackCharge &charge, const GlobalVector &momemtum, const GlobalPoint &position, double &xc, double &yc, double &r, const MagneticField &magField) const
void TangentApproachInRPhi::circleParameters ( const TrackCharge charge,
const GlobalVector momemtum,
const GlobalPoint position,
double &  xc,
double &  yc,
double &  r,
const MagneticField magField 
) const
private

temporary code, to be replaced by call to curvature() when bug is fixed.

end of temporary code

Definition at line 165 of file TangentApproachInRPhi.cc.

References funct::abs(), funct::cos(), MagneticField::inTesla(), PV3DBase< T, PVType, FrameType >::phi(), phi(), funct::sin(), PV3DBase< T, PVType, FrameType >::transverse(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

171 {
172 
173  // compute radius of circle
177 // double bz = MagneticField::inInverseGeV(position).z();
178  double bz = magField.inTesla(position).z() * 2.99792458e-3;
179 
180  // signed_r directed towards circle center, along F_Lorentz = q*v X B
181  double signed_r = charge*momentum.transverse() / bz;
182  r = abs(signed_r);
186  // compute centre of circle
187  double phi = momentum.phi();
188  xc = signed_r*sin(phi) + position.x();
189  yc = -signed_r*cos(phi) + position.y();
190 
191 }
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T y() const
Definition: PV3DBase.h:63
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T x() const
Definition: PV3DBase.h:62
virtual TangentApproachInRPhi* TangentApproachInRPhi::clone ( void  ) const
inlinevirtual

Clone method

Implements ClosestApproachOnHelices.

Definition at line 51 of file TangentApproachInRPhi.h.

References TangentApproachInRPhi().

51  {
52  return new TangentApproachInRPhi(* this);
53  }
GlobalPoint TangentApproachInRPhi::crossingPoint ( ) const
virtual

arithmetic mean of the two points of closest approach

Implements ClosestApproachOnHelices.

Definition at line 47 of file TangentApproachInRPhi.cc.

References Exception.

Referenced by ConversionProducer::preselectTrackPair().

48 {
49  if (!status_)
50  throw cms::Exception("TrackingTools/PatternTools","TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
51  return GlobalPoint((posA.x() + posB.x())/2.,
52  (posA.y() + posB.y())/2.,
53  (posA.z() + posB.z())/2.);
54 }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
T z() const
Definition: PV3DBase.h:64
T x() const
Definition: PV3DBase.h:62
float TangentApproachInRPhi::distance ( ) const
virtual

distance between the two points of closest approach in 3D

Implements ClosestApproachOnHelices.

Definition at line 57 of file TangentApproachInRPhi.cc.

References Exception.

58 {
59  if (!status_)
60  throw cms::Exception("TrackingTools/PatternTools","TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
61  return (posB - posA).mag();
62 }
float TangentApproachInRPhi::perpdist ( ) const

signed distance between two points of closest approach in r-phi plane (-ive if circles intersect)

Definition at line 64 of file TangentApproachInRPhi.cc.

References Exception, and perp().

Referenced by ConversionProducer::preselectTrackPair().

65 {
66  if (!status_)
67  throw cms::Exception("TrackingTools/PatternTools","TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
68 
69  float perpdist = (posB - posA).perp();
70 
71  if (intersection_) {
72  perpdist = -perpdist;
73  }
74 
75  return perpdist;
76 
77 }
T perp() const
Magnitude of transverse component.
pair< GlobalPoint, GlobalPoint > TangentApproachInRPhi::points ( ) const
virtual

Returns the two PCA on the trajectories.

Implements ClosestApproachOnHelices.

Definition at line 38 of file TangentApproachInRPhi.cc.

References Exception.

39 {
40  if (!status_)
41  throw cms::Exception("TrackingTools/PatternTools","TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
42  return pair<GlobalPoint, GlobalPoint> (posA, posB);
43 }
virtual bool TangentApproachInRPhi::status ( void  ) const
inlinevirtual

Implements ClosestApproachOnHelices.

Definition at line 27 of file TangentApproachInRPhi.h.

References status_.

Referenced by ConversionProducer::preselectTrackPair().

pair< GlobalTrajectoryParameters, GlobalTrajectoryParameters > TangentApproachInRPhi::trajectoryParameters ( ) const

Returns not only the points, but the full GlobalTrajectoryParemeters at the points of closest approach

Definition at line 126 of file TangentApproachInRPhi.cc.

References Exception, and run_regression::ret.

Referenced by ConversionProducer::preselectTrackPair().

127 {
128  if (!status_)
129  throw cms::Exception("TrackingTools/PatternTools","TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
130  pair <GlobalTrajectoryParameters, GlobalTrajectoryParameters>
133  return ret;
134 }
std::pair< GlobalTrajectoryParameters, GlobalTrajectoryParameters > trajectoryParameters() const
GlobalTrajectoryParameters paramA
GlobalTrajectoryParameters paramB
GlobalTrajectoryParameters TangentApproachInRPhi::trajectoryParameters ( const GlobalPoint newpt,
const GlobalTrajectoryParameters oldpar 
) const
private

Definition at line 136 of file TangentApproachInRPhi.cc.

References GlobalTrajectoryParameters::charge(), GlobalTrajectoryParameters::magneticField(), GlobalTrajectoryParameters::momentum(), GlobalTrajectoryParameters::position(), alignCSCRings::r, mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

138 {
139  // First we need the centers of the circles.
140  double xc, yc, r;
141  circleParameters( oldgtp.charge(), oldgtp.momentum(),
142  oldgtp.position(), xc, yc, r, oldgtp.magneticField() );
143 
144  // now we do a translation, move the center of circle to (0,0,0).
145  double dx1 = oldgtp.position().x() - xc;
146  double dy1 = oldgtp.position().y() - yc;
147  double dx2 = newpt.x() - xc;
148  double dy2 = newpt.y() - yc;
149 
150  // now for the angles:
151  double cosphi = ( dx1 * dx2 + dy1 * dy2 ) /
152  ( sqrt ( dx1 * dx1 + dy1 * dy1 ) * sqrt ( dx2 * dx2 + dy2 * dy2 ));
153  double sinphi = - oldgtp.charge() * sqrt ( 1 - cosphi * cosphi );
154 
155  // Finally, the new momenta:
156  double px = cosphi * oldgtp.momentum().x() - sinphi * oldgtp.momentum().y();
157  double py = sinphi * oldgtp.momentum().x() + cosphi * oldgtp.momentum().y();
158 
159  GlobalVector vta ( px, py, oldgtp.momentum().z() );
160  GlobalTrajectoryParameters gta( newpt , vta , oldgtp.charge(), &(oldgtp.magneticField()) );
161  return gta;
162 }
T y() const
Definition: PV3DBase.h:63
T sqrt(T t)
Definition: SSEVec.h:48
void circleParameters(const TrackCharge &charge, const GlobalVector &momemtum, const GlobalPoint &position, double &xc, double &yc, double &r, const MagneticField &magField) const
T x() const
Definition: PV3DBase.h:62
int TangentApproachInRPhi::transverseCoord ( double  cxa,
double  cya,
double  ra,
double  cxb,
double  cyb,
double  rb,
double &  xg1,
double &  yg1,
double &  xg2,
double &  yg2 
) const
private

Definition at line 195 of file TangentApproachInRPhi.cc.

References funct::abs(), jetcorrextractor::sign(), mathSSE::sqrt(), and findQualityFiles::v.

199 {
200  int flag = 0;
201  double x1, y1, x2, y2;
202 
203  // new reference frame with origin in (cxa, cya) and x-axis
204  // directed from (cxa, cya) to (cxb, cyb)
205 
206  double d_ab = sqrt((cxb - cxa)*(cxb - cxa) + (cyb - cya)*(cyb - cya));
207  if (d_ab == 0) { // concentric circles
208  return 0;
209  }
210  // elements of rotation matrix
211  double u = (cxb - cxa) / d_ab;
212  double v = (cyb - cya) / d_ab;
213 
214  // conditions for circle intersection
215  if (d_ab <= ra + rb && d_ab >= abs(rb - ra)) {
216 
217  // circles cross each other
218 // flag = 1;
219 //
220 // // triangle (ra, rb, d_ab)
221 // double cosphi = (ra*ra - rb*rb + d_ab*d_ab) / (2*ra*d_ab);
222 // double sinphi2 = 1. - cosphi*cosphi;
223 // if (sinphi2 < 0.) { sinphi2 = 0.; cosphi = 1.; }
224 //
225 // // intersection points in new frame
226 // double sinphi = sqrt(sinphi2);
227 // x1 = ra*cosphi; y1 = ra*sinphi; x2 = x1; y2 = -y1;
228 
229  //circles cross each other, but take tangent points anyway
230  flag = 1;
231 
232  // points of closest approach in new frame
233  // are on line between 2 centers
234  x1 = ra; y1 = 0; x2 = d_ab - rb; y2 = 0;
235 
236 
237  }
238  else if (d_ab > ra + rb) {
239 
240  // circles are external to each other
241  flag = 2;
242 
243  // points of closest approach in new frame
244  // are on line between 2 centers
245  x1 = ra; y1 = 0; x2 = d_ab - rb; y2 = 0;
246  }
247  else if (d_ab < abs(rb - ra)) {
248 
249  // circles are inside each other
250  flag = 2;
251 
252  // points of closest approach in new frame are on line between 2 centers
253  // choose 2 closest points
254  double sign = 1.;
255  if (ra <= rb) sign = -1.;
256  x1 = sign*ra; y1 = 0; x2 = d_ab + sign*rb; y2 = 0;
257  }
258  else {
259  return 0;
260  }
261 
262  // intersection points in global frame, transverse plane
263  xg1 = u*x1 - v*y1 + cxa; yg1 = v*x1 + u*y1 + cya;
264  xg2 = u*x2 - v*y2 + cxa; yg2 = v*x2 + u*y2 + cya;
265 
266  return flag;
267 }
double sign(double x)
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double TangentApproachInRPhi::zCoord ( const GlobalVector mom,
const GlobalPoint pos,
double  r,
double  xc,
double  yc,
double  xg,
double  yg 
) const
private

Definition at line 271 of file TangentApproachInRPhi.cc.

References funct::abs(), phi(), alignCSCRings::r, mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::transverse(), PV3DBase< T, PVType, FrameType >::x(), x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

275 {
276 
277  // starting point
278  double x = pos.x(); double y = pos.y(); double z = pos.z();
279 
280  double px = mom.x(); double py = mom.y(); double pz = mom.z();
281 
282  // rotation angle phi from starting point to crossing point (absolute value)
283  // -- compute sin(phi/2) if phi smaller than pi/4,
284  // -- cos(phi) if phi larger than pi/4
285  double phi = 0.;
286  double sinHalfPhi = sqrt((x-xg)*(x-xg) + (y-yg)*(y-yg))/(2*r);
287  if (sinHalfPhi < 0.383) { // sin(pi/8)
288  phi = 2*asin(sinHalfPhi);
289  }
290  else {
291  double cosPhi = ((x-xc)*(xg-xc) + (y-yc)*(yg-yc))/(r*r);
292  if (std::abs(cosPhi) > 1) cosPhi = (cosPhi > 0 ? 1 : -1);
293  phi = abs(acos(cosPhi));
294  }
295  // -- sign of phi
296  double signPhi = ((x - xc)*(yg - yc) - (xg - xc)*(y - yc) > 0) ? 1. : -1.;
297 
298  // sign of track angular momentum
299  // if rotation is along angular momentum, delta z is along pz
300  double signOmega = ((x - xc)*py - (y - yc)*px > 0) ? 1. : -1.;
301 
302  // delta z
303  // -- |dz| = |cos(theta) * path along helix|
304  // = |cos(theta) * arc length along circle / sin(theta)|
305  double dz = signPhi*signOmega*(pz/mom.transverse())*phi*r;
306 
307  return z + dz;
308 }
T y() const
Definition: PV3DBase.h:63
T transverse() const
Definition: PV3DBase.h:73
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T x() const
Definition: PV3DBase.h:62

Member Data Documentation

bool TangentApproachInRPhi::intersection_
private

Definition at line 100 of file TangentApproachInRPhi.h.

Referenced by TangentApproachInRPhi().

GlobalTrajectoryParameters TangentApproachInRPhi::paramA
private

Definition at line 99 of file TangentApproachInRPhi.h.

GlobalTrajectoryParameters TangentApproachInRPhi::paramB
private

Definition at line 99 of file TangentApproachInRPhi.h.

GlobalPoint TangentApproachInRPhi::posA
private

Definition at line 98 of file TangentApproachInRPhi.h.

GlobalPoint TangentApproachInRPhi::posB
private

Definition at line 98 of file TangentApproachInRPhi.h.

bool TangentApproachInRPhi::status_
private

Definition at line 97 of file TangentApproachInRPhi.h.

Referenced by status(), and TangentApproachInRPhi().