CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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

bool calculate (const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
 
bool calculate (const FreeTrajectoryState &sta, const FreeTrajectoryState &stb) override
 
TangentApproachInRPhiclone () const override
 
GlobalPoint crossingPoint () const override
 
float distance () const override
 
float perpdist () const
 
std::pair< GlobalPoint,
GlobalPoint
points () const override
 
bool status () const override
 
 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 17 of file TangentApproachInRPhi.h.

References intersection_, and status_.

Referenced by clone().

17  {
18  status_ = false;
19  intersection_ = false;
20  }

Member Function Documentation

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

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().

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

Implements ClosestApproachOnHelices.

Definition at line 21 of file TangentApproachInRPhi.cc.

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

21  {
22  TrackCharge chargeA = sta.charge();
23  TrackCharge chargeB = stb.charge();
24  GlobalVector momentumA = sta.momentum();
25  GlobalVector momentumB = stb.momentum();
26  GlobalPoint positionA = sta.position();
27  GlobalPoint positionB = stb.position();
28  paramA = sta.parameters();
29  paramB = stb.parameters();
30 
31  return calculate(chargeA, momentumA, positionA, chargeB, momentumB, positionB, sta.parameters().magneticField());
32 }
const GlobalTrajectoryParameters & parameters() const
GlobalTrajectoryParameters paramA
TrackCharge charge() const
int TrackCharge
Definition: TrackCharge.h:4
GlobalTrajectoryParameters paramB
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
GlobalVector momentum() const
GlobalPoint position() const
const MagneticField & magneticField() const
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 73 of file TangentApproachInRPhi.cc.

79  {
80  // centres and radii of track circles
81  double xca, yca, ra;
82  circleParameters(chargeA, momentumA, positionA, xca, yca, ra, magField);
83  double xcb, ycb, rb;
84  circleParameters(chargeB, momentumB, positionB, xcb, ycb, rb, magField);
85 
86  // points of closest approach in transverse plane
87  double xg1, yg1, xg2, yg2;
88  int flag = transverseCoord(xca, yca, ra, xcb, ycb, rb, xg1, yg1, xg2, yg2);
89  if (flag == 0) {
90  status_ = false;
91  return false;
92  }
93 
94  double xga, yga, zga, xgb, ygb, zgb;
95 
96  if (flag == 1) {
97  intersection_ = true;
98  } else {
99  intersection_ = false;
100  }
101 
102  // one point of closest approach on each track in transverse plane
103  xga = xg1;
104  yga = yg1;
105  zga = zCoord(momentumA, positionA, ra, xca, yca, xga, yga);
106  xgb = xg2;
107  ygb = yg2;
108  zgb = zCoord(momentumB, positionB, rb, xcb, ycb, xgb, ygb);
109 
110  posA = GlobalPoint(xga, yga, zga);
111  posB = GlobalPoint(xgb, ygb, zgb);
112  status_ = true;
113  return true;
114 }
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 151 of file TangentApproachInRPhi.cc.

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

157  {
158  // compute radius of circle
162  // double bz = MagneticField::inInverseGeV(position).z();
163  double bz = magField.inTesla(position).z() * 2.99792458e-3;
164 
165  // signed_r directed towards circle center, along F_Lorentz = q*v X B
166  double signed_r = charge * momentum.transverse() / bz;
167  r = abs(signed_r);
171  // compute centre of circle
172  double phi = momentum.phi();
173  xc = signed_r * sin(phi) + position.x();
174  yc = -signed_r * cos(phi) + position.y();
175 }
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:60
T z() const
Definition: PV3DBase.h:61
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:59
TangentApproachInRPhi* TangentApproachInRPhi::clone ( void  ) const
inlineoverridevirtual

Clone method

Implements ClosestApproachOnHelices.

Definition at line 49 of file TangentApproachInRPhi.h.

References TangentApproachInRPhi().

49 { return new TangentApproachInRPhi(*this); }
GlobalPoint TangentApproachInRPhi::crossingPoint ( ) const
overridevirtual

arithmetic mean of the two points of closest approach

Implements ClosestApproachOnHelices.

Definition at line 42 of file TangentApproachInRPhi.cc.

References Exception.

Referenced by ConversionProducer::preselectTrackPair().

42  {
43  if (!status_)
44  throw cms::Exception(
45  "TrackingTools/PatternTools",
46  "TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
47  return GlobalPoint((posA.x() + posB.x()) / 2., (posA.y() + posB.y()) / 2., (posA.z() + posB.z()) / 2.);
48 }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:60
T z() const
Definition: PV3DBase.h:61
T x() const
Definition: PV3DBase.h:59
float TangentApproachInRPhi::distance ( ) const
overridevirtual

distance between the two points of closest approach in 3D

Implements ClosestApproachOnHelices.

Definition at line 50 of file TangentApproachInRPhi.cc.

References Exception.

50  {
51  if (!status_)
52  throw cms::Exception(
53  "TrackingTools/PatternTools",
54  "TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
55  return (posB - posA).mag();
56 }
float TangentApproachInRPhi::perpdist ( ) const

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

Definition at line 58 of file TangentApproachInRPhi.cc.

References Exception, and perp().

Referenced by ConversionProducer::preselectTrackPair().

58  {
59  if (!status_)
60  throw cms::Exception(
61  "TrackingTools/PatternTools",
62  "TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
63 
64  float perpdist = (posB - posA).perp();
65 
66  if (intersection_) {
67  perpdist = -perpdist;
68  }
69 
70  return perpdist;
71 }
T perp() const
Magnitude of transverse component.
pair< GlobalPoint, GlobalPoint > TangentApproachInRPhi::points ( ) const
overridevirtual

Returns the two PCA on the trajectories.

Implements ClosestApproachOnHelices.

Definition at line 34 of file TangentApproachInRPhi.cc.

References Exception.

34  {
35  if (!status_)
36  throw cms::Exception(
37  "TrackingTools/PatternTools",
38  "TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
39  return pair<GlobalPoint, GlobalPoint>(posA, posB);
40 }
bool TangentApproachInRPhi::status ( void  ) const
inlineoverridevirtual

Implements ClosestApproachOnHelices.

Definition at line 26 of file TangentApproachInRPhi.h.

References status_.

Referenced by ConversionProducer::preselectTrackPair().

26 { return status_; }
pair< GlobalTrajectoryParameters, GlobalTrajectoryParameters > TangentApproachInRPhi::trajectoryParameters ( ) const

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

Definition at line 116 of file TangentApproachInRPhi.cc.

References Exception, and runTheMatrix::ret.

Referenced by ConversionProducer::preselectTrackPair().

116  {
117  if (!status_)
118  throw cms::Exception(
119  "TrackingTools/PatternTools",
120  "TangentApproachInRPhi::could not compute track crossing. Check status before calling this method!");
121  pair<GlobalTrajectoryParameters, GlobalTrajectoryParameters> ret(trajectoryParameters(posA, paramA),
123  return ret;
124 }
tuple ret
prodAgent to be discontinued
std::pair< GlobalTrajectoryParameters, GlobalTrajectoryParameters > trajectoryParameters() const
GlobalTrajectoryParameters paramA
GlobalTrajectoryParameters paramB
GlobalTrajectoryParameters TangentApproachInRPhi::trajectoryParameters ( const GlobalPoint newpt,
const GlobalTrajectoryParameters oldpar 
) const
private

Definition at line 126 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().

127  {
128  // First we need the centers of the circles.
129  double xc, yc, r;
130  circleParameters(oldgtp.charge(), oldgtp.momentum(), oldgtp.position(), xc, yc, r, oldgtp.magneticField());
131 
132  // now we do a translation, move the center of circle to (0,0,0).
133  double dx1 = oldgtp.position().x() - xc;
134  double dy1 = oldgtp.position().y() - yc;
135  double dx2 = newpt.x() - xc;
136  double dy2 = newpt.y() - yc;
137 
138  // now for the angles:
139  double cosphi = (dx1 * dx2 + dy1 * dy2) / (sqrt(dx1 * dx1 + dy1 * dy1) * sqrt(dx2 * dx2 + dy2 * dy2));
140  double sinphi = -oldgtp.charge() * sqrt(1 - cosphi * cosphi);
141 
142  // Finally, the new momenta:
143  double px = cosphi * oldgtp.momentum().x() - sinphi * oldgtp.momentum().y();
144  double py = sinphi * oldgtp.momentum().x() + cosphi * oldgtp.momentum().y();
145 
146  GlobalVector vta(px, py, oldgtp.momentum().z());
147  GlobalTrajectoryParameters gta(newpt, vta, oldgtp.charge(), &(oldgtp.magneticField()));
148  return gta;
149 }
T y() const
Definition: PV3DBase.h:60
T sqrt(T t)
Definition: SSEVec.h:19
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:59
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 177 of file TangentApproachInRPhi.cc.

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

186  {
187  int flag = 0;
188  double x1, y1, x2, y2;
189 
190  // new reference frame with origin in (cxa, cya) and x-axis
191  // directed from (cxa, cya) to (cxb, cyb)
192 
193  double d_ab = sqrt((cxb - cxa) * (cxb - cxa) + (cyb - cya) * (cyb - cya));
194  if (d_ab == 0) { // concentric circles
195  return 0;
196  }
197  // elements of rotation matrix
198  double u = (cxb - cxa) / d_ab;
199  double v = (cyb - cya) / d_ab;
200 
201  // conditions for circle intersection
202  if (d_ab <= ra + rb && d_ab >= abs(rb - ra)) {
203  // circles cross each other
204  // flag = 1;
205  //
206  // // triangle (ra, rb, d_ab)
207  // double cosphi = (ra*ra - rb*rb + d_ab*d_ab) / (2*ra*d_ab);
208  // double sinphi2 = 1. - cosphi*cosphi;
209  // if (sinphi2 < 0.) { sinphi2 = 0.; cosphi = 1.; }
210  //
211  // // intersection points in new frame
212  // double sinphi = sqrt(sinphi2);
213  // x1 = ra*cosphi; y1 = ra*sinphi; x2 = x1; y2 = -y1;
214 
215  //circles cross each other, but take tangent points anyway
216  flag = 1;
217 
218  // points of closest approach in new frame
219  // are on line between 2 centers
220  x1 = ra;
221  y1 = 0;
222  x2 = d_ab - rb;
223  y2 = 0;
224 
225  } else if (d_ab > ra + rb) {
226  // circles are external to each other
227  flag = 2;
228 
229  // points of closest approach in new frame
230  // are on line between 2 centers
231  x1 = ra;
232  y1 = 0;
233  x2 = d_ab - rb;
234  y2 = 0;
235  } else if (d_ab < abs(rb - ra)) {
236  // circles are inside each other
237  flag = 2;
238 
239  // points of closest approach in new frame are on line between 2 centers
240  // choose 2 closest points
241  double sign = 1.;
242  if (ra <= rb)
243  sign = -1.;
244  x1 = sign * ra;
245  y1 = 0;
246  x2 = d_ab + sign * rb;
247  y2 = 0;
248  } else {
249  return 0;
250  }
251 
252  // intersection points in global frame, transverse plane
253  xg1 = u * x1 - v * y1 + cxa;
254  yg1 = v * x1 + u * y1 + cya;
255  xg2 = u * x2 - v * y2 + cxa;
256  yg2 = v * x2 + u * y2 + cya;
257 
258  return flag;
259 }
double sign(double x)
T sqrt(T t)
Definition: SSEVec.h:19
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 261 of file TangentApproachInRPhi.cc.

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

262  {
263  // starting point
264  double x = pos.x();
265  double y = pos.y();
266  double z = pos.z();
267 
268  double px = mom.x();
269  double py = mom.y();
270  double pz = mom.z();
271 
272  // rotation angle phi from starting point to crossing point (absolute value)
273  // -- compute sin(phi/2) if phi smaller than pi/4,
274  // -- cos(phi) if phi larger than pi/4
275  double phi = 0.;
276  double sinHalfPhi = sqrt((x - xg) * (x - xg) + (y - yg) * (y - yg)) / (2 * r);
277  if (sinHalfPhi < 0.383) { // sin(pi/8)
278  phi = 2 * asin(sinHalfPhi);
279  } else {
280  double cosPhi = ((x - xc) * (xg - xc) + (y - yc) * (yg - yc)) / (r * r);
281  if (std::abs(cosPhi) > 1)
282  cosPhi = (cosPhi > 0 ? 1 : -1);
283  phi = abs(acos(cosPhi));
284  }
285  // -- sign of phi
286  double signPhi = ((x - xc) * (yg - yc) - (xg - xc) * (y - yc) > 0) ? 1. : -1.;
287 
288  // sign of track angular momentum
289  // if rotation is along angular momentum, delta z is along pz
290  double signOmega = ((x - xc) * py - (y - yc) * px > 0) ? 1. : -1.;
291 
292  // delta z
293  // -- |dz| = |cos(theta) * path along helix|
294  // = |cos(theta) * arc length along circle / sin(theta)|
295  double dz = signPhi * signOmega * (pz / mom.transverse()) * phi * r;
296 
297  return z + dz;
298 }
T y() const
Definition: PV3DBase.h:60
T transverse() const
Definition: PV3DBase.h:70
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T x() const
Definition: PV3DBase.h:59

Member Data Documentation

bool TangentApproachInRPhi::intersection_
private

Definition at line 102 of file TangentApproachInRPhi.h.

Referenced by TangentApproachInRPhi().

GlobalTrajectoryParameters TangentApproachInRPhi::paramA
private

Definition at line 101 of file TangentApproachInRPhi.h.

GlobalTrajectoryParameters TangentApproachInRPhi::paramB
private

Definition at line 101 of file TangentApproachInRPhi.h.

GlobalPoint TangentApproachInRPhi::posA
private

Definition at line 100 of file TangentApproachInRPhi.h.

GlobalPoint TangentApproachInRPhi::posB
private

Definition at line 100 of file TangentApproachInRPhi.h.

bool TangentApproachInRPhi::status_
private

Definition at line 99 of file TangentApproachInRPhi.h.

Referenced by status(), and TangentApproachInRPhi().