StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StPhysicalHelix.cc
1 /***************************************************************************
2  *
3  * $Id: StPhysicalHelix.cc,v 1.8 2007/08/20 23:24:54 perev Exp $
4  *
5  * Author: Brian Lasiuk, Sep 1997
6  ***************************************************************************
7  *
8  * Description:
9  * Parametrization of a physical helix.
10  *
11  ***************************************************************************
12  *
13  * $Log: StPhysicalHelix.cc,v $
14  * Revision 1.8 2007/08/20 23:24:54 perev
15  * BugFix #1027
16  *
17  * Revision 1.7 2005/07/06 18:49:56 fisyak
18  * Replace StHelixD, StLorentzVectorD,StLorentzVectorF,StMatrixD,StMatrixF,StPhysicalHelixD,StThreeVectorD,StThreeVectorF by templated version
19  *
20  * Revision 1.6 2003/09/02 17:59:35 perev
21  * gcc 3.2 updates + WarnOff
22  *
23  * Revision 1.5 2002/06/21 17:49:26 genevb
24  * Some minor speed improvements
25  *
26  * Revision 1.4 2002/02/20 00:56:23 ullrich
27  * Added methods to calculate signed DCA.
28  *
29  * Revision 1.3 1999/02/24 11:42:18 ullrich
30  * Fixed bug in momentum().
31  *
32  * Revision 1.2 1999/02/12 01:01:04 wenaus
33  * Fix bug in momentum calculation
34  *
35  * Revision 1.1 1999/01/30 03:59:04 fisyak
36  * Root Version of StarClassLibrary
37  *
38  * Revision 1.1 1999/01/23 00:29:21 ullrich
39  * Initial Revision
40  *
41  **************************************************************************/
42 #include <math.h>
43 #include "StHelix.hh"
44 #include "StPhysicalHelix.hh"
45 #include "PhysicalConstants.h"
46 #include "SystemOfUnits.h"
47 #ifdef __ROOT__
48 ClassImpT(StPhysicalHelix,double);
49 #endif
50 StPhysicalHelix::StPhysicalHelix(){}
51 
52 StPhysicalHelix::~StPhysicalHelix() { /* nop */ }
53 
54 StPhysicalHelix::StPhysicalHelix(const StThreeVector<double>& p,
55  const StThreeVector<double>& o,
56  double B, double q)
57 {
58  mH = (q*B <= 0) ? 1 : -1;
59  if(p.y() == 0 && p.x() == 0)
60  setPhase((M_PI/4)*(1-2.*mH));
61  else
62  setPhase(atan2(p.y(),p.x())-mH*M_PI/2);
63  setDipAngle(atan2(p.z(),p.perp()));
64  mOrigin = o;
65 
66 #ifndef ST_NO_NAMESPACES
67  {
68  using namespace units;
69 #endif
70  setCurvature(fabs((c_light*nanosecond/meter*q*B/tesla)/
71  (abs(p)/GeV*mCosDipAngle)/meter));
72 #ifndef ST_NO_NAMESPACES
73  }
74 #endif
75 }
76 
77 StPhysicalHelix::StPhysicalHelix(double c, double d, double phase,
78  const StThreeVector<double>& o, int h)
79  : StHelix(c, d, phase, o, h) { /* nop */}
80 
81 
82 StThreeVector<double> StPhysicalHelix::momentum(double B) const
83 {
84  if (mSingularity)
85  return(StThreeVector<double>(0,0,0));
86  else {
87 #ifndef ST_NO_NAMESPACES
88  {
89  using namespace units;
90 #endif
91  double pt = GeV*fabs(c_light*nanosecond/meter*B/tesla)/(fabs(mCurvature)*meter);
92 
93  return (StThreeVector<double>(pt*cos(mPhase+mH*M_PI/2), // pos part pos field
94  pt*sin(mPhase+mH*M_PI/2),
95  pt*tan(mDipAngle)));
96 #ifndef ST_NO_NAMESPACES
97  }
98 #endif
99  }
100 }
101 
102 StThreeVector<double> StPhysicalHelix::momentumAt(double S, double B) const
103 {
104  // Obtain phase-shifted momentum from phase-shift of origin
105  StPhysicalHelix tmp(*this);
106  tmp.moveOrigin(S);
107  return tmp.momentum(B);
108 }
109 
110 int StPhysicalHelix::charge(double B) const
111 {
112  return (B > 0 ? -mH : mH);
113 }
114 
115 double StPhysicalHelix::geometricSignedDistance(double x, double y)
116 {
117  // Geometric signed distance
118  double thePath = this->pathLength(x,y);
119  StThreeVector<double> DCA2dPosition = this->at(thePath);
120  DCA2dPosition.setZ(0);
121  StThreeVector<double> position(x,y,0);
122  StThreeVector<double> DCAVec = (DCA2dPosition-position);
123  StThreeVector<double> momVec;
124  // Deal with straight tracks
125  if (this->mSingularity) {
126  momVec = this->at(1)- this->at(0);
127  momVec.setZ(0);
128  }
129  else {
130  momVec = this->momentumAt(thePath,1./tesla); // Don't care about Bmag. Helicity is what matters.
131  momVec.setZ(0);
132  }
133 
134  double cross = DCAVec.x()*momVec.y() - DCAVec.y()*momVec.x();
135  double theSign = (cross>=0) ? 1. : -1.;
136  return theSign*DCAVec.perp();
137 }
138 
139 double StPhysicalHelix::curvatureSignedDistance(double x, double y)
140 {
141  // Protect against mH = 0 or zero field
142  if (this->mSingularity || abs(this->mH)<=0) {
143  return (this->geometricSignedDistance(x,y));
144  }
145  else {
146  return (this->geometricSignedDistance(x,y))/(this->mH);
147  }
148 
149 }
150 
151 double StPhysicalHelix::geometricSignedDistance(const StThreeVector<double>& pos)
152 {
153  double sdca2d = this->geometricSignedDistance(pos.x(),pos.y());
154  double theSign = (sdca2d>=0) ? 1. : -1.;
155  return (this->distance(pos))*theSign;
156 }
157 
158 double StPhysicalHelix::curvatureSignedDistance(const StThreeVector<double>& pos)
159 {
160  double sdca2d = this->curvatureSignedDistance(pos.x(),pos.y());
161  double theSign = (sdca2d>=0) ? 1. : -1.;
162  return (this->distance(pos))*theSign;
163 }
int h() const
y-center of circle in xy-plane
Definition: StHelix.hh:174
double x(double s) const
coordinates of helix at point s
Definition: StHelix.hh:182
double distance(const StThreeVector< double > &p, bool scanPeriods=true) const
minimal distance between point and helix
Definition: StHelix.cc:240
pair< double, double > pathLength(double r) const
path length at given r (cylindrical r)
Definition: StHelix.cc:351
void setPhase(double)
performs also various checks
Definition: StHelix.cc:191
double phase() const
1/R in xy-plane
Definition: StHelix.hh:180