CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
fastsim::HelixTrajectory Class Reference

Mathematical representation of a helix. More...

#include <HelixTrajectory.h>

Inheritance diagram for fastsim::HelixTrajectory:
fastsim::Trajectory

Public Member Functions

bool crosses (const BarrelSimplifiedGeometry &layer) const override
 Check if an intersection of the trajectory with a barrel layer exists. More...
 
double getRadParticle (double phi) const
 Return distance of particle from center of the detector if it was at given angle phi of the helix. More...
 
 HelixTrajectory (const Particle &particle, double magneticFieldZ)
 Constructor. More...
 
void move (double deltaTimeC) override
 Move the particle along the helix trajectory for a given time. More...
 
double nextCrossingTimeC (const BarrelSimplifiedGeometry &layer, bool onLayer=false) const override
 Return delta time (t*c) of the next intersection of trajectory and barrel layer. More...
 
- Public Member Functions inherited from fastsim::Trajectory
const math::XYZTLorentzVectorgetMomentum ()
 Simple getter: return momentum of the particle that was used to create trajectory. More...
 
const math::XYZTLorentzVectorgetPosition ()
 Simple getter: return position of the particle that was used to create trajectory. More...
 
double nextCrossingTimeC (const SimplifiedGeometry &layer, bool onLayer=false) const
 Return delta time (t*c) of the next intersection of trajectory and generic layer. More...
 
double nextCrossingTimeC (const ForwardSimplifiedGeometry &layer, bool onLayer=false) const
 Return delta time (t*c) of the next intersection of trajectory and forward layer. More...
 
virtual ~Trajectory ()
 

Private Attributes

const double centerR_
 Distance of the center of the helix from the center of the tracker. More...
 
const double centerX_
 X-coordinate of the center of the helix. More...
 
const double centerY_
 Y-coordinate of the center of the helix. More...
 
const double maxR_
 The maximum distance of the helix from the center of the tracker. More...
 
const double minR_
 The minimal distance of the helix from the center of the tracker. More...
 
const double phi_
 
const double phiSpeed_
 The angular speed of the particle on the helix trajectory. More...
 
const double radius_
 The radius of the helix. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from fastsim::Trajectory
static std::unique_ptr< TrajectorycreateTrajectory (const fastsim::Particle &particle, const double magneticFieldZ)
 Calls constructor of derived classes. More...
 
- Protected Member Functions inherited from fastsim::Trajectory
 Trajectory (const fastsim::Particle &particle)
 Constructor of base class. More...
 
- Protected Attributes inherited from fastsim::Trajectory
math::XYZTLorentzVector momentum_
 momentum of the particle that was used to create trajectory More...
 
math::XYZTLorentzVector position_
 position of the particle that was used to create trajectory More...
 

Detailed Description

Mathematical representation of a helix.

Reflects trajectory of a charged particle in a magnetic field. The trajectory is defined by cylindrical coordinates (see definition of variables for more information).

Definition at line 21 of file HelixTrajectory.h.

Constructor & Destructor Documentation

fastsim::HelixTrajectory::HelixTrajectory ( const Particle particle,
double  magneticFieldZ 
)

Constructor.

The magnetic field is (to good approximation) constant between two tracker layers (and only in Z-direction).

Parameters
particleA (usually charged) particle.
magneticFieldZThe magnetic field.

Definition at line 15 of file HelixTrajectory.cc.

16  : Trajectory(particle)
17  // exact: r = gamma*beta*m_0*c / (q*e*B) = p_T / (q * e * B)
18  // momentum in units of GeV/c: r = p_T * 10^9 / (c * q * B)
19  // in cmssw units: r = p_T / (c * 10^-4 * q * B)
20  , radius_(std::abs(momentum_.Pt() / (fastsim::Constants::speedOfLight * 1e-4 * particle.charge() * magneticFieldZ)))
21  , phi_(std::atan(momentum_.Py()/momentum_.Px()) + (momentum_.Px()*particle.charge() < 0 ? 3.*M_PI/2. : M_PI/2. ))
22  // maybe consider (for -pi/2<x<pi/2)
23  // cos(atan(x)) = 1 / sqrt(x^2+1)
24  // -> cos(atan(x) + pi/2) = - x / sqrt(x^2+1)
25  // -> cos(atan(x) +3*pi/2) = + x / sqrt(x^2+1)
26  // sin(atan(x)) = x / sqrt(x^2+1)
27  // -> sin(atan(x) + pi/2) = + 1 / sqrt(x^2+1)
28  // -> sin(atan(x) +3*pi/2) = - 1 / sqrt(x^2+1)
29  , centerX_(position_.X() - radius_ * (momentum_.Py()/momentum_.Px()) / std::sqrt((momentum_.Py()/momentum_.Px())*(momentum_.Py()/momentum_.Px())+1) * (momentum_.Px()*particle.charge() < 0 ? 1. : -1.))
30  , centerY_(position_.Y() - radius_ * 1 / std::sqrt((momentum_.Py()/momentum_.Px())*(momentum_.Py()/momentum_.Px())+1) * (momentum_.Px()*particle.charge() < 0 ? -1. : 1.))
31  //, centerX_(position_.X() - radius_*std::cos(phi_))
32  //, centerY_(position_.Y() - radius_*std::sin(phi_))
36  // omega = q * e * B / (gamma * m) = q * e *B / (E / c^2) = q * e * B * c^2 / E
37  // omega: negative for negative q -> seems to be what we want.
38  // energy in units of GeV: omega = q * B * c^2 / (E * 10^9)
39  // in cmssw units: omega[1/ns] = q * B * c^2 * 10^-4 / E
40  , phiSpeed_(-particle.charge() * magneticFieldZ * fastsim::Constants::speedOfLight * fastsim::Constants::speedOfLight * 1e-4 / momentum_.E())
41 {;}
const double minR_
The minimal distance of the helix from the center of the tracker.
const double radius_
The radius of the helix.
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:16
const double centerY_
Y-coordinate of the center of the helix.
math::XYZTLorentzVector position_
position of the particle that was used to create trajectory
Definition: Trajectory.h:95
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Trajectory(const fastsim::Particle &particle)
Constructor of base class.
Definition: Trajectory.cc:12
const double maxR_
The maximum distance of the helix from the center of the tracker.
const double phiSpeed_
The angular speed of the particle on the helix trajectory.
#define M_PI
const double centerX_
X-coordinate of the center of the helix.
const double centerR_
Distance of the center of the helix from the center of the tracker.
math::XYZTLorentzVector momentum_
momentum of the particle that was used to create trajectory
Definition: Trajectory.h:96

Member Function Documentation

bool fastsim::HelixTrajectory::crosses ( const BarrelSimplifiedGeometry layer) const
overridevirtual

Check if an intersection of the trajectory with a barrel layer exists.

Parameters
layerA barrel layer.

Implements fastsim::Trajectory.

Definition at line 43 of file HelixTrajectory.cc.

References fastsim::BarrelSimplifiedGeometry::getRadius(), maxR_, and minR_.

Referenced by nextCrossingTimeC().

44 {
45  return (minR_ < layer.getRadius() && maxR_ > layer.getRadius());
46 }
const double minR_
The minimal distance of the helix from the center of the tracker.
const double maxR_
The maximum distance of the helix from the center of the tracker.
double fastsim::HelixTrajectory::getRadParticle ( double  phi) const

Return distance of particle from center of the detector if it was at given angle phi of the helix.

Parameters
phiangle of the helix

Definition at line 210 of file HelixTrajectory.cc.

References centerX_, centerY_, funct::cos(), radius_, funct::sin(), and mathSSE::sqrt().

Referenced by nextCrossingTimeC().

211 {
212  return sqrt((centerX_ + radius_*std::cos(phi))*(centerX_ + radius_*std::cos(phi)) + (centerY_ + radius_*std::sin(phi))*(centerY_ + radius_*std::sin(phi)));
213 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const double radius_
The radius of the helix.
const double centerY_
Y-coordinate of the center of the helix.
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const double centerX_
X-coordinate of the center of the helix.
void fastsim::HelixTrajectory::move ( double  deltaTimeC)
overridevirtual

Move the particle along the helix trajectory for a given time.

Parameters
deltaTimeCTime in units of t*c..

Implements fastsim::Trajectory.

Definition at line 191 of file HelixTrajectory.cc.

References centerX_, centerY_, funct::cos(), hiPixelPairStep_cff::deltaPhi, fastsim::Trajectory::momentum_, phi_, phiSpeed_, fastsim::Trajectory::position_, radius_, funct::sin(), and fastsim::Constants::speedOfLight.

Referenced by Vispa.Gui.WidgetContainer.WidgetContainer::autosize(), Vispa.Gui.VispaWidget.VispaWidget::dragWidget(), Vispa.Gui.VispaWidget.VispaWidget::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

192 {
193  double deltaT = deltaTimeC/fastsim::Constants::speedOfLight;
194  double deltaPhi = phiSpeed_*deltaT;
195  position_.SetXYZT(
196  centerX_ + radius_*std::cos(phi_ + deltaPhi),
197  centerY_ + radius_*std::sin(phi_ + deltaPhi),
198  position_.Z() + momentum_.Z()/momentum_.E()*deltaTimeC,
199  position_.T() + deltaT);
200  // Rotation defined by
201  // x' = x cos θ - y sin θ
202  // y' = x sin θ + y cos θ
203  momentum_.SetXYZT(
204  momentum_.X()*std::cos(deltaPhi) - momentum_.Y()*std::sin(deltaPhi),
205  momentum_.X()*std::sin(deltaPhi) + momentum_.Y()*std::cos(deltaPhi),
206  momentum_.Z(),
207  momentum_.E());
208 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const double radius_
The radius of the helix.
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:16
const double centerY_
Y-coordinate of the center of the helix.
math::XYZTLorentzVector position_
position of the particle that was used to create trajectory
Definition: Trajectory.h:95
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const double phiSpeed_
The angular speed of the particle on the helix trajectory.
const double centerX_
X-coordinate of the center of the helix.
math::XYZTLorentzVector momentum_
momentum of the particle that was used to create trajectory
Definition: Trajectory.h:96
double fastsim::HelixTrajectory::nextCrossingTimeC ( const BarrelSimplifiedGeometry layer,
bool  onLayer = false 
) const
overridevirtual

Return delta time (t*c) of the next intersection of trajectory and barrel layer.

This function solves the quadratic equation (basically intersection of two circles with a given radius) in order to calculate the moment in time when the particle's trajectory intersects with a given barrel layer.

Parameters
layerA barrel layer.
onLayerSpecify if the particle already is on the layer (in that case the second solution has to be picked).
Returns
t*c [ns * cm/ns] of next intersection (-1 if there is none).

Implements fastsim::Trajectory.

Definition at line 48 of file HelixTrajectory.cc.

References a, funct::abs(), b, EnergyCorrector::c, centerX_, centerY_, funct::cos(), crosses(), delta, MillePedeFileConverter_cfg::e, F(), callgraph::G, fastsim::BarrelSimplifiedGeometry::getRadius(), getRadParticle(), M_PI, min(), fastsim::Trajectory::momentum_, fastsim::StraightTrajectory::nextCrossingTimeC(), phi_, phiSpeed_, fastsim::Trajectory::position_, radius_, funct::sin(), fastsim::Constants::speedOfLight, and mathSSE::sqrt().

49 {
50  if(!crosses(layer)) return -1;
51 
52  // solve the following equation for sin(phi)
53  // (x^2 + y^2 = R_L^2) (1) the layer
54  // x = x_c + R_H*cos(phi) (2) the helix in the xy plane
55  // y = y_c + R_H*sin(phi) (3) the helix in the xy plane
56  // with
57  // R_L: the radius of the layer
58  // x_c,y_c the center of the helix in xy plane
59  // R_H, the radius of the helix
60  // phi, the phase of the helix
61  //
62  // substitute (2) and (3) in (1)
63  // =>
64  // x_c^2 + 2*x_c*R_H*cos(phi) + R_H^2*cos^2(phi)
65  // + y_c^2 + 2*y_c*R_H*sin(phi) + R_H^2*sin^2(phi)
66  // = R_L^2
67  // =>
68  // (x_c^2 + y_c^2 + R_H^2 - R_L^2) + (2*y_c*R_H)*sin(phi) = -(2*x_c*R_H)*cos(phi)
69  //
70  // rewrite
71  // E + F *sin(phi) = G *cos(phi)
72  // =>
73  // E^2 + 2*E*F*sin(phi) + F^2*sin^2(phi) = G^2*(1-sin^2(phi))
74  // rearrange
75  // sin^2(phi)*(F^2 + G^2) + sin(phi)*(2*E*F) + (E^2 - G^2) = 0
76  //
77  // rewrite
78  // sin^2(phi)* a + sin(phi)* b + c = 0
79  // => sin(phi) = (-b +/- sqrt(b^2 - 4*ac)) / (2*a)
80  // with
81  // a = F^2 + G^2
82  // b = 2*E*F
83  // c = E^2 - G^2
84 
85  double E = centerX_*centerX_ + centerY_*centerY_ + radius_*radius_ - layer.getRadius()*layer.getRadius();
86  double F = 2*centerY_*radius_;
87  double G = 2*centerX_*radius_;
88 
89  double a = F*F + G*G;
90  double b = 2*E*F;
91  double c = E*E - G*G;
92 
93  double delta = b*b - 4*a*c;
94 
95  // case of no solution
96  if(delta < 0)
97  {
98  // Should not be reached: Full Propagation does always have a solution "if(crosses(layer)) == -1"
99  // Even if particle is outside all layers -> can turn around in magnetic field
100  return -1;
101  }
102 
103  // Uses a numerically more stable procedure:
104  // https://people.csail.mit.edu/bkph/articles/Quadratics.pdf
105  double sqrtDelta = sqrt(delta);
106  double phi1 = 0, phi2 = 0;
107  if(b < 0){
108  phi1 = std::asin((2.*c) / (-b + sqrtDelta));
109  phi2 = std::asin((-b + sqrtDelta) / (2.*a));
110  }else{
111  phi1 = std::asin((-b - sqrtDelta) / (2.*a));
112  phi2 = std::asin((2.*c) / (-b - sqrtDelta));
113  }
114 
115  // asin is ambiguous, make sure to have the right solution
116  if(std::abs(layer.getRadius() - getRadParticle(phi1)) > 1.0e-2){
117  phi1 = - phi1 + M_PI;
118  }
119  if(std::abs(layer.getRadius() - getRadParticle(phi2)) > 1.0e-2){
120  phi2 = - phi2 + M_PI;
121  }
122 
123  // another ambiguity
124  if(phi1 < 0){
125  phi1 += 2. * M_PI;
126  }
127  if(phi2 < 0){
128  phi2 += 2. * M_PI;
129  }
130 
131  // find the corresponding times when the intersection occurs
132  // make sure they are positive
133  double t1 = (phi1 - phi_)/phiSpeed_;
134  while(t1 < 0)
135  {
136  t1 += 2*M_PI/std::abs(phiSpeed_);
137  }
138  double t2 = (phi2 - phi_)/phiSpeed_;
139  while(t2 < 0)
140  {
141  t2 += 2*M_PI/std::abs(phiSpeed_);
142  }
143 
144  // Check if propagation successful (numerical reasons): both solutions (phi1, phi2) have to be on the layer (same radius)
145  // Can happen due to numerical instabilities of geometrical function (if momentum is almost parallel to x/y axis)
146  // Get crossingTimeC from StraightTrajectory as good approximation
147  if(std::abs(layer.getRadius() - getRadParticle(phi1)) > 1.0e-2
148  || std::abs(layer.getRadius() - getRadParticle(phi2)) > 1.0e-2)
149  {
150  StraightTrajectory traj(*this);
151  return traj.nextCrossingTimeC(layer, onLayer);
152  }
153 
154  // if the particle is already on the layer, we need to make sure the 2nd solution is picked.
155  // happens if particle turns around in the magnetic field instead of hitting the next layer
156  if(onLayer)
157  {
158  bool particleMovesInwards = momentum_.X()*position_.X() + momentum_.Y()*position_.Y() < 0;
159 
160  double posX1 = centerX_ + radius_*std::cos(phi1);
161  double posY1 = centerY_ + radius_*std::sin(phi1);
162  double momX1 = momentum_.X()*std::cos(phi1 - phi_) - momentum_.Y()*std::sin(phi1 - phi_);
163  double momY1 = momentum_.X()*std::sin(phi1 - phi_) + momentum_.Y()*std::cos(phi1 - phi_);
164  bool particleMovesInwards1 = momX1*posX1 + momY1*posY1 < 0;
165 
166  double posX2 = centerX_ + radius_*std::cos(phi2);
167  double posY2 = centerY_ + radius_*std::sin(phi2);
168  double momX2 = momentum_.X()*std::cos(phi2 - phi_) - momentum_.Y()*std::sin(phi2 - phi_);
169  double momY2 = momentum_.X()*std::sin(phi2 - phi_) + momentum_.Y()*std::cos(phi2 - phi_);
170  bool particleMovesInwards2 = momX2*posX2 + momY2*posY2 < 0;
171 
172  if(particleMovesInwards1 != particleMovesInwards)
173  {
175  }
176  else if(particleMovesInwards2 != particleMovesInwards)
177  {
179  }
180  // try to catch numerical issues again..
181  else
182  {
183  return -1;
184  }
185  }
186 
188 }
dbl * delta
Definition: mlp_gen.cc:36
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const double radius_
The radius of the helix.
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:16
const double centerY_
Y-coordinate of the center of the helix.
double getRadParticle(double phi) const
Return distance of particle from center of the detector if it was at given angle phi of the helix...
math::XYZTLorentzVector position_
position of the particle that was used to create trajectory
Definition: Trajectory.h:95
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
const double phiSpeed_
The angular speed of the particle on the helix trajectory.
#define M_PI
bool crosses(const BarrelSimplifiedGeometry &layer) const override
Check if an intersection of the trajectory with a barrel layer exists.
const double centerX_
X-coordinate of the center of the helix.
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
math::XYZTLorentzVector momentum_
momentum of the particle that was used to create trajectory
Definition: Trajectory.h:96
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281

Member Data Documentation

const double fastsim::HelixTrajectory::centerR_
private

Distance of the center of the helix from the center of the tracker.

Definition at line 66 of file HelixTrajectory.h.

const double fastsim::HelixTrajectory::centerX_
private

X-coordinate of the center of the helix.

Definition at line 64 of file HelixTrajectory.h.

Referenced by getRadParticle(), move(), and nextCrossingTimeC().

const double fastsim::HelixTrajectory::centerY_
private

Y-coordinate of the center of the helix.

Definition at line 65 of file HelixTrajectory.h.

Referenced by getRadParticle(), move(), and nextCrossingTimeC().

const double fastsim::HelixTrajectory::maxR_
private

The maximum distance of the helix from the center of the tracker.

Definition at line 68 of file HelixTrajectory.h.

Referenced by crosses().

const double fastsim::HelixTrajectory::minR_
private

The minimal distance of the helix from the center of the tracker.

Definition at line 67 of file HelixTrajectory.h.

Referenced by crosses().

const double fastsim::HelixTrajectory::phi_
private

The angle of the particle alone the helix. Ranges from 0 to 2PI: 0 corresponds to the positive X direction, phi increases counterclockwise

Definition at line 62 of file HelixTrajectory.h.

Referenced by move(), and nextCrossingTimeC().

const double fastsim::HelixTrajectory::phiSpeed_
private

The angular speed of the particle on the helix trajectory.

Definition at line 69 of file HelixTrajectory.h.

Referenced by move(), and nextCrossingTimeC().

const double fastsim::HelixTrajectory::radius_
private

The radius of the helix.

Definition at line 61 of file HelixTrajectory.h.

Referenced by getRadParticle(), move(), and nextCrossingTimeC().