Public Member Functions | Private Attributes | Static Private Attributes | List of all members
trkf::SurfYZLine Class Reference

#include <SurfYZLine.h>

Inheritance diagram for trkf::SurfYZLine:
trkf::SurfLine trkf::Surface trkf::SurfWireLine

Public Member Functions

 SurfYZLine ()
 Default constructor. More...
 
 SurfYZLine (double x0, double y0, double z0, double phi)
 Initializing constructor. More...
 
virtual ~SurfYZLine ()
 Destructor. More...
 
double x0 () const
 X origin. More...
 
double y0 () const
 Y origin. More...
 
double z0 () const
 Z origin. More...
 
double phi () const
 Rotation angle about x-axis. More...
 
virtual Surfaceclone () const
 Clone method. More...
 
virtual bool isTrackValid (const TrackVector &vec) const
 Surface-specific tests of validity of track parameters. More...
 
virtual void toLocal (const double xyz[3], double uvw[3]) const
 Transform global to local coordinates. More...
 
virtual void toGlobal (const double uvw[3], double xyz[3]) const
 Transform local to global coordinates. More...
 
virtual TrackVector getDiff (const TrackVector &vec1, const TrackVector &vec2) const
 Calculate difference of two track parameter vectors. More...
 
virtual void getPosition (const TrackVector &vec, double xyz[3]) const
 Get position of track. More...
 
virtual void getMomentum (const TrackVector &vec, double mom[3], TrackDirection dir=UNKNOWN) const
 Get momentum vector of track. More...
 
virtual bool isParallel (const Surface &surf) const
 Test whether two surfaces are parallel, within tolerance. More...
 
virtual double distanceTo (const Surface &surf) const
 Find perpendicular distance to a parallel surface. More...
 
virtual bool isEqual (const Surface &surf) const
 Test two surfaces for equality, within tolerance. More...
 
virtual std::ostream & Print (std::ostream &out) const
 Printout. More...
 
- Public Member Functions inherited from trkf::SurfLine
 SurfLine ()
 Default constructor. More...
 
virtual ~SurfLine ()
 Destructor. More...
 
double PointingError (const TrackVector &vec, const TrackError &err) const
 Get pointing error of track. More...
 
void getStartingError (TrackError &err) const
 Get starting error matrix for Kalman filter. More...
 
- Public Member Functions inherited from trkf::Surface
 Surface ()
 Default constructor. More...
 
virtual ~Surface ()
 Destructor. More...
 
virtual TrackDirection getDirection (const TrackVector &, TrackDirection dir=UNKNOWN) const
 Get direction of track (default UNKNOWN). More...
 

Private Attributes

double fX0
 X origin. More...
 
double fY0
 Y origin. More...
 
double fZ0
 Z origin. More...
 
double fPhi
 Rotation angle about x-axis. More...
 

Static Private Attributes

static double fPhiTolerance = 1.e-10
 Phi tolerance for parallel. More...
 
static double fSepTolerance = 1.e-6
 Separation tolerance for equal. More...
 

Additional Inherited Members

- Public Types inherited from trkf::Surface
enum  TrackDirection { FORWARD, BACKWARD, UNKNOWN }
 Track direction enum. More...
 

Detailed Description

Definition at line 78 of file SurfYZLine.h.

Constructor & Destructor Documentation

trkf::SurfYZLine::SurfYZLine ( )

Default constructor.

Definition at line 25 of file SurfYZLine.cxx.

25  :
26  fX0(0.),
27  fY0(0.),
28  fZ0(0.),
29  fPhi(0.)
30  {}
double fX0
X origin.
Definition: SurfYZLine.h:140
double fZ0
Z origin.
Definition: SurfYZLine.h:142
double fY0
Y origin.
Definition: SurfYZLine.h:141
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
trkf::SurfYZLine::SurfYZLine ( double  x0,
double  y0,
double  z0,
double  phi 
)

Initializing constructor.

Initializing constructor.

Arguments:

x0, y0, z0 - Global coordinates of local origin. phi - Rotation angle about x-axis.

Definition at line 39 of file SurfYZLine.cxx.

39  :
40  fX0(x0),
41  fY0(y0),
42  fZ0(z0),
43  fPhi(phi)
44  {}
double fX0
X origin.
Definition: SurfYZLine.h:140
double z0() const
Z origin.
Definition: SurfYZLine.h:94
double x0() const
X origin.
Definition: SurfYZLine.h:92
double fZ0
Z origin.
Definition: SurfYZLine.h:142
double fY0
Y origin.
Definition: SurfYZLine.h:141
double phi() const
Rotation angle about x-axis.
Definition: SurfYZLine.h:95
double y0() const
Y origin.
Definition: SurfYZLine.h:93
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
trkf::SurfYZLine::~SurfYZLine ( )
virtual

Destructor.

Definition at line 47 of file SurfYZLine.cxx.

48  {}

Member Function Documentation

Surface * trkf::SurfYZLine::clone ( ) const
virtual

Clone method.

Implements trkf::Surface.

Definition at line 51 of file SurfYZLine.cxx.

52  {
53  return new SurfYZLine(*this);
54  }
SurfYZLine()
Default constructor.
Definition: SurfYZLine.cxx:25
double trkf::SurfYZLine::distanceTo ( const Surface surf) const
virtual

Find perpendicular distance to a parallel surface.

Find perpendicular distance to a parallel surface.

Throw an exception if the other surface is not parallel.

Arguments:

surf - Other surface.

Returned value: Distance.

Implements trkf::Surface.

Definition at line 237 of file SurfYZLine.cxx.

238  {
239  // Check if the other surface is parallel to this one.
240 
241  bool parallel = isParallel(surf);
242  if(!parallel)
243  throw cet::exception("SurfYZLine") << "Attempt to find distance to non-parallel surface.\n";
244 
245  // Find the origin of the other surface in global coordinates,
246  // then convert to our local coordinates.
247 
248  double otheruvw[3] = {0., 0., 0.};
249  double xyz[3];
250  double myuvw[3];
251  surf.toGlobal(otheruvw, xyz);
252  toLocal(xyz, myuvw);
253 
254  // Distance of v-axis to other surface origin.
255 
256  return std::hypot(myuvw[0], myuvw[2]);
257  }
virtual void toLocal(const double xyz[3], double uvw[3]) const
Transform global to local coordinates.
Definition: SurfYZLine.cxx:71
std::enable_if_t< std::is_arithmetic_v< T >, T > hypot(T x, T y)
Definition: hypot.h:60
virtual bool isParallel(const Surface &surf) const
Test whether two surfaces are parallel, within tolerance.
Definition: SurfYZLine.cxx:208
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
TrackVector trkf::SurfYZLine::getDiff ( const TrackVector vec1,
const TrackVector vec2 
) const
virtual

Calculate difference of two track parameter vectors.

Calculate difference of two track parameter vectors, taking into account phi wrap.

Arguments:

vec1 - First vector. vec2 - Second vector.

Returns: vec1 - vec2

Reimplemented from trkf::Surface.

Definition at line 117 of file SurfYZLine.cxx.

118  {
119  TrackVector result = vec1 - vec2;
120  while(result(2) <= -TMath::Pi())
121  result(2) += TMath::TwoPi();
122  while(result(2) > TMath::Pi())
123  result(2) -= TMath::TwoPi();
124  return result;
125  }
static QCString result
KVector< 5 >::type TrackVector
Track state vector, dimension 5.
void trkf::SurfYZLine::getMomentum ( const TrackVector vec,
double  mom[3],
TrackDirection  dir = UNKNOWN 
) const
virtual

Get momentum vector of track.

Get momentum vector of track.

Arguments:

vec - Track state vector. mom - Momentum vector in global coordinate system. dir - Track direction (ignored).

Implements trkf::Surface.

Definition at line 161 of file SurfYZLine.cxx.

163  {
164 
165  // Get momentum.
166 
167  double invp = std::abs(vec(4));
168  double p = 1. / std::max(invp, 1.e-3); // Capped at 1000. GeV/c.
169 
170  // Get track direction parameters.
171 
172  double phi = vec(2);
173  double eta = vec(3);
174 
175  double sinphi = std::sin(phi);
176  double cosphi = std::cos(phi);
177  double sh = 1./std::cosh(eta); // sech(eta)
178  double th = std::tanh(eta);
179 
180  // Calculate momentum vector in local coordinate system.
181 
182  double pu = p * cosphi * sh;
183  double pv = p * th;
184  double pw = p * sinphi * sh;
185 
186  // Rotate momentum to global coordinte system.
187 
188  double sinfphi = std::sin(fPhi);
189  double cosfphi = std::cos(fPhi);
190 
191  mom[0] = pu;
192  mom[1] = pv * cosfphi - pw * sinfphi;
193  mom[2] = pv * sinfphi + pw * cosfphi;
194 
195  return;
196  }
T abs(T value)
const double e
p
Definition: test.py:223
static int max(int a, int b)
double phi() const
Rotation angle about x-axis.
Definition: SurfYZLine.h:95
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
void trkf::SurfYZLine::getPosition ( const TrackVector vec,
double  xyz[3] 
) const
virtual

Get position of track.

Get position of track.

Arguments:

vec - Track state vector. xyz - Position in global coordinate system.

Implements trkf::Surface.

Definition at line 134 of file SurfYZLine.cxx.

135  {
136  // Get position in local coordinate system.
137 
138  double phi = vec(2);
139  double sinphi = std::sin(phi);
140  double cosphi = std::cos(phi);
141 
142  double uvw[3];
143  uvw[0] = -vec(0) * sinphi;
144  uvw[1] = vec(1);
145  uvw[2] = vec(0) * cosphi;
146 
147  // Transform to global coordinate system.
148 
149  toGlobal(uvw, xyz);
150  return;
151  }
virtual void toGlobal(const double uvw[3], double xyz[3]) const
Transform local to global coordinates.
Definition: SurfYZLine.cxx:93
double phi() const
Rotation angle about x-axis.
Definition: SurfYZLine.h:95
bool trkf::SurfYZLine::isEqual ( const Surface surf) const
virtual

Test two surfaces for equality, within tolerance.

Test two surfaces for equality, within tolerance. Here equal is defined as having all surface parameters the same, not just having the surfaces coincide spatially, so that the local coordinate systems are the same between the two surfaces.

Arguments:

surf - Other surface.

Returned values: true if equal.

Implements trkf::Surface.

Definition at line 270 of file SurfYZLine.cxx.

271  {
272  bool result = false;
273 
274  // Test if the other surface is a SurfYZLine.
275 
276  const SurfYZLine* psurf = dynamic_cast<const SurfYZLine*>(&surf);
277  if(psurf != 0) {
278 
279  // Test whether surface parameters are the same within tolerance.
280 
281  double delta_phi = TVector2::Phi_mpi_pi(fPhi - psurf->phi());
282  double dx = fX0 - psurf->x0();
283  double dy = fY0 - psurf->y0();
284  double dz = fZ0 - psurf->z0();
285  if(std::abs(delta_phi) <= fPhiTolerance &&
286  std::abs(dx) <= fSepTolerance &&
287  std::abs(dy) <= fSepTolerance &&
288  std::abs(dz) <= fSepTolerance)
289  result = true;
290  }
291  return result;
292  }
static QCString result
T abs(T value)
double fX0
X origin.
Definition: SurfYZLine.h:140
static double fPhiTolerance
Phi tolerance for parallel.
Definition: SurfYZLine.h:135
static double fSepTolerance
Separation tolerance for equal.
Definition: SurfYZLine.h:136
double fZ0
Z origin.
Definition: SurfYZLine.h:142
double fY0
Y origin.
Definition: SurfYZLine.h:141
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
SurfYZLine()
Default constructor.
Definition: SurfYZLine.cxx:25
bool trkf::SurfYZLine::isParallel ( const Surface surf) const
virtual

Test whether two surfaces are parallel, within tolerance.

Test whether two surfaces are parallel, within tolerance. This method will only return true if the other surface is a SurfYZLine.

Arguments:

surf - Other surface.

Returned value: true if parallel.

Implements trkf::Surface.

Definition at line 208 of file SurfYZLine.cxx.

209  {
210  bool result = false;
211 
212  // Test if the other surface is a SurfYZLine.
213 
214  const SurfYZLine* psurf = dynamic_cast<const SurfYZLine*>(&surf);
215  if(psurf != 0) {
216 
217  // Test whether surface angle parameters are the same
218  // within tolerance.
219 
220  double delta_phi = TVector2::Phi_mpi_pi(fPhi - psurf->phi());
221  if(std::abs(delta_phi) <= fPhiTolerance)
222  result = true;
223  }
224  return result;
225  }
static QCString result
T abs(T value)
static double fPhiTolerance
Phi tolerance for parallel.
Definition: SurfYZLine.h:135
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
SurfYZLine()
Default constructor.
Definition: SurfYZLine.cxx:25
bool trkf::SurfYZLine::isTrackValid ( const TrackVector vec) const
virtual

Surface-specific tests of validity of track parameters.

Implements trkf::Surface.

Definition at line 57 of file SurfYZLine.cxx.

58  {
59  // Limit allowed range of eta parameter.
60 
61  return std::abs(vec(3)) < 10.;
62  }
T abs(T value)
double trkf::SurfYZLine::phi ( ) const
inline

Rotation angle about x-axis.

Definition at line 95 of file SurfYZLine.h.

std::ostream & trkf::SurfYZLine::Print ( std::ostream &  out) const
virtual

Printout.

Implements trkf::Surface.

Definition at line 295 of file SurfYZLine.cxx.

296  {
297  out << "SurfYZLine{ x0=" << fX0 << ", y0=" << fY0 << ", z0=" << fZ0 << ", phi=" << fPhi << "}";
298  return out;
299  }
double fX0
X origin.
Definition: SurfYZLine.h:140
double fZ0
Z origin.
Definition: SurfYZLine.h:142
double fY0
Y origin.
Definition: SurfYZLine.h:141
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
void trkf::SurfYZLine::toGlobal ( const double  uvw[3],
double  xyz[3] 
) const
virtual

Transform local to global coordinates.

Transform local to global coordinates.

Arguments:

uvw - Cartesian coordinates in local coordinate system. xyz - Cartesian coordinates in global coordinate system.

Implements trkf::Surface.

Definition at line 93 of file SurfYZLine.cxx.

94  {
95  double sinphi = std::sin(fPhi);
96  double cosphi = std::cos(fPhi);
97 
98  // x = x0 + u
99  xyz[0] = fX0 + uvw[0];
100 
101  // y = y0 + v*cos(phi) - w*sin(phi)
102  xyz[1] = fY0 + uvw[1] * cosphi - uvw[2] * sinphi;
103 
104  // z = z0 + v*sin(phi) + w*cos(phi)
105  xyz[2] = fZ0 + uvw[1] * sinphi + uvw[2] * cosphi;
106  }
double fX0
X origin.
Definition: SurfYZLine.h:140
double fZ0
Z origin.
Definition: SurfYZLine.h:142
double fY0
Y origin.
Definition: SurfYZLine.h:141
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
void trkf::SurfYZLine::toLocal ( const double  xyz[3],
double  uvw[3] 
) const
virtual

Transform global to local coordinates.

Transform global to local coordinates.

Arguments:

xyz - Cartesian coordinates in global coordinate system. uvw - Cartesian coordinates in local coordinate system.

Implements trkf::Surface.

Definition at line 71 of file SurfYZLine.cxx.

72  {
73  double sinphi = std::sin(fPhi);
74  double cosphi = std::cos(fPhi);
75 
76  // u = x-x0
77  uvw[0] = xyz[0] - fX0;
78 
79  // v = (y-y0)*cos(phi) + (z-z0)*sin(phi)
80  uvw[1] = (xyz[1] - fY0) * cosphi + (xyz[2] - fZ0) * sinphi;
81 
82  // w = -(y-y0)*sin(phi) + (z-z0)*cos(phi)
83  uvw[2] = -(xyz[1] - fY0) * sinphi + (xyz[2] - fZ0) * cosphi;
84  }
double fX0
X origin.
Definition: SurfYZLine.h:140
double fZ0
Z origin.
Definition: SurfYZLine.h:142
double fY0
Y origin.
Definition: SurfYZLine.h:141
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:143
double trkf::SurfYZLine::x0 ( ) const
inline

X origin.

Definition at line 92 of file SurfYZLine.h.

double trkf::SurfYZLine::y0 ( ) const
inline

Y origin.

Definition at line 93 of file SurfYZLine.h.

double trkf::SurfYZLine::z0 ( ) const
inline

Z origin.

Definition at line 94 of file SurfYZLine.h.

Member Data Documentation

double trkf::SurfYZLine::fPhi
private

Rotation angle about x-axis.

Definition at line 143 of file SurfYZLine.h.

double trkf::SurfYZLine::fPhiTolerance = 1.e-10
staticprivate

Phi tolerance for parallel.

Definition at line 135 of file SurfYZLine.h.

double trkf::SurfYZLine::fSepTolerance = 1.e-6
staticprivate

Separation tolerance for equal.

Definition at line 136 of file SurfYZLine.h.

double trkf::SurfYZLine::fX0
private

X origin.

Definition at line 140 of file SurfYZLine.h.

double trkf::SurfYZLine::fY0
private

Y origin.

Definition at line 141 of file SurfYZLine.h.

double trkf::SurfYZLine::fZ0
private

Z origin.

Definition at line 142 of file SurfYZLine.h.


The documentation for this class was generated from the following files: