LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
trkf::SurfYZLine Class Reference

#include "SurfYZLine.h"

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

Public Types

enum  TrackDirection { FORWARD, BACKWARD, UNKNOWN }
 Track direction enum. More...
 

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

Detailed Description

Definition at line 78 of file SurfYZLine.h.

Member Enumeration Documentation

Track direction enum.

Enumerator
FORWARD 
BACKWARD 
UNKNOWN 

Definition at line 54 of file Surface.h.

Constructor & Destructor Documentation

trkf::SurfYZLine::SurfYZLine ( )

Default constructor.

Definition at line 25 of file SurfYZLine.cxx.

Referenced by clone(), and trkf::SurfWireLine::SurfWireLine().

25 : fX0(0.), fY0(0.), fZ0(0.), fPhi(0.) {}
double fX0
X origin.
Definition: SurfYZLine.h:138
double fZ0
Z origin.
Definition: SurfYZLine.h:140
double fY0
Y origin.
Definition: SurfYZLine.h:139
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
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 34 of file SurfYZLine.cxx.

35  : fX0(x0), fY0(y0), fZ0(z0), fPhi(phi)
36  {}
double fX0
X origin.
Definition: SurfYZLine.h:138
double z0() const
Z origin.
Definition: SurfYZLine.h:92
double x0() const
X origin.
Definition: SurfYZLine.h:90
double fZ0
Z origin.
Definition: SurfYZLine.h:140
double fY0
Y origin.
Definition: SurfYZLine.h:139
double phi() const
Rotation angle about x-axis.
Definition: SurfYZLine.h:93
double y0() const
Y origin.
Definition: SurfYZLine.h:91
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
trkf::SurfYZLine::~SurfYZLine ( )
virtual

Destructor.

Definition at line 39 of file SurfYZLine.cxx.

39 {}

Member Function Documentation

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

Clone method.

Implements trkf::Surface.

Definition at line 42 of file SurfYZLine.cxx.

References SurfYZLine().

Referenced by phi().

43  {
44  return new SurfYZLine(*this);
45  }
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 226 of file SurfYZLine.cxx.

References isParallel(), trkf::Surface::toGlobal(), and toLocal().

Referenced by phi().

227  {
228  // Check if the other surface is parallel to this one.
229 
230  bool parallel = isParallel(surf);
231  if (!parallel)
232  throw cet::exception("SurfYZLine") << "Attempt to find distance to non-parallel surface.\n";
233 
234  // Find the origin of the other surface in global coordinates,
235  // then convert to our local coordinates.
236 
237  double otheruvw[3] = {0., 0., 0.};
238  double xyz[3];
239  double myuvw[3];
240  surf.toGlobal(otheruvw, xyz);
241  toLocal(xyz, myuvw);
242 
243  // Distance of v-axis to other surface origin.
244 
245  return std::hypot(myuvw[0], myuvw[2]);
246  }
virtual void toLocal(const double xyz[3], double uvw[3]) const
Transform global to local coordinates.
Definition: SurfYZLine.cxx:62
virtual bool isParallel(const Surface &surf) const
Test whether two surfaces are parallel, within tolerance.
Definition: SurfYZLine.cxx:198
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 108 of file SurfYZLine.cxx.

Referenced by phi().

109  {
110  TrackVector result = vec1 - vec2;
111  while (result(2) <= -TMath::Pi())
112  result(2) += TMath::TwoPi();
113  while (result(2) > TMath::Pi())
114  result(2) -= TMath::TwoPi();
115  return result;
116  }
KVector< 5 >::type TrackVector
Track state vector, dimension 5.
virtual TrackDirection trkf::Surface::getDirection ( const TrackVector ,
TrackDirection  dir = UNKNOWN 
) const
inlinevirtualinherited
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 152 of file SurfYZLine.cxx.

References util::abs(), e, fPhi, and phi().

Referenced by phi().

153  {
154 
155  // Get momentum.
156 
157  double invp = std::abs(vec(4));
158  double p = 1. / std::max(invp, 1.e-3); // Capped at 1000. GeV/c.
159 
160  // Get track direction parameters.
161 
162  double phi = vec(2);
163  double eta = vec(3);
164 
165  double sinphi = std::sin(phi);
166  double cosphi = std::cos(phi);
167  double sh = 1. / std::cosh(eta); // sech(eta)
168  double th = std::tanh(eta);
169 
170  // Calculate momentum vector in local coordinate system.
171 
172  double pu = p * cosphi * sh;
173  double pv = p * th;
174  double pw = p * sinphi * sh;
175 
176  // Rotate momentum to global coordinte system.
177 
178  double sinfphi = std::sin(fPhi);
179  double cosfphi = std::cos(fPhi);
180 
181  mom[0] = pu;
182  mom[1] = pv * cosfphi - pw * sinfphi;
183  mom[2] = pv * sinfphi + pw * cosfphi;
184 
185  return;
186  }
constexpr auto abs(T v)
Returns the absolute value of the argument.
double phi() const
Rotation angle about x-axis.
Definition: SurfYZLine.h:93
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
Float_t e
Definition: plot.C:35
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 125 of file SurfYZLine.cxx.

References phi(), and toGlobal().

Referenced by phi().

126  {
127  // Get position in local coordinate system.
128 
129  double phi = vec(2);
130  double sinphi = std::sin(phi);
131  double cosphi = std::cos(phi);
132 
133  double uvw[3];
134  uvw[0] = -vec(0) * sinphi;
135  uvw[1] = vec(1);
136  uvw[2] = vec(0) * cosphi;
137 
138  // Transform to global coordinate system.
139 
140  toGlobal(uvw, xyz);
141  return;
142  }
virtual void toGlobal(const double uvw[3], double xyz[3]) const
Transform local to global coordinates.
Definition: SurfYZLine.cxx:84
double phi() const
Rotation angle about x-axis.
Definition: SurfYZLine.h:93
void trkf::SurfLine::getStartingError ( TrackError err) const
virtualinherited

Get starting error matrix for Kalman filter.

Get starting error matrix for Kalman filter.

Arguments:

err - Error matrix.

Implements trkf::Surface.

Definition at line 90 of file SurfLine.cxx.

91  {
92  err.resize(5, false);
93  err.clear();
94  err(0, 0) = 1000.;
95  err(1, 1) = 1000.;
96  err(2, 2) = 10.;
97  err(3, 3) = 10.;
98  err(4, 4) = 10.;
99  }
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 259 of file SurfYZLine.cxx.

References util::abs(), fPhi, fPhiTolerance, fSepTolerance, fX0, fY0, fZ0, phi(), x0(), y0(), and z0().

Referenced by trkf::KHitWireLine::KHitWireLine(), and phi().

260  {
261  bool result = false;
262 
263  // Test if the other surface is a SurfYZLine.
264 
265  const SurfYZLine* psurf = dynamic_cast<const SurfYZLine*>(&surf);
266  if (psurf != 0) {
267 
268  // Test whether surface parameters are the same within tolerance.
269 
270  double delta_phi = TVector2::Phi_mpi_pi(fPhi - psurf->phi());
271  double dx = fX0 - psurf->x0();
272  double dy = fY0 - psurf->y0();
273  double dz = fZ0 - psurf->z0();
274  if (std::abs(delta_phi) <= fPhiTolerance && std::abs(dx) <= fSepTolerance &&
275  std::abs(dy) <= fSepTolerance && std::abs(dz) <= fSepTolerance)
276  result = true;
277  }
278  return result;
279  }
constexpr auto abs(T v)
Returns the absolute value of the argument.
double fX0
X origin.
Definition: SurfYZLine.h:138
static double fPhiTolerance
Phi tolerance for parallel.
Definition: SurfYZLine.h:133
static double fSepTolerance
Separation tolerance for equal.
Definition: SurfYZLine.h:134
double fZ0
Z origin.
Definition: SurfYZLine.h:140
double fY0
Y origin.
Definition: SurfYZLine.h:139
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
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 198 of file SurfYZLine.cxx.

References util::abs(), fPhi, fPhiTolerance, and phi().

Referenced by distanceTo(), and phi().

199  {
200  bool result = false;
201 
202  // Test if the other surface is a SurfYZLine.
203 
204  const SurfYZLine* psurf = dynamic_cast<const SurfYZLine*>(&surf);
205  if (psurf != 0) {
206 
207  // Test whether surface angle parameters are the same
208  // within tolerance.
209 
210  double delta_phi = TVector2::Phi_mpi_pi(fPhi - psurf->phi());
211  if (std::abs(delta_phi) <= fPhiTolerance) result = true;
212  }
213  return result;
214  }
constexpr auto abs(T v)
Returns the absolute value of the argument.
static double fPhiTolerance
Phi tolerance for parallel.
Definition: SurfYZLine.h:133
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
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 48 of file SurfYZLine.cxx.

References util::abs().

Referenced by phi().

49  {
50  // Limit allowed range of eta parameter.
51 
52  return std::abs(vec(3)) < 10.;
53  }
constexpr auto abs(T v)
Returns the absolute value of the argument.
double trkf::SurfLine::PointingError ( const TrackVector vec,
const TrackError err 
) const
virtualinherited

Get pointing error of track.

Get pointing error of track.

Arguments:

vec - Track parameters. err - Track error matrix.

Returns: Pointing error.

This method calculates the track pointing error based on the slope track paramers and errors (parameters 2 and 3).

Implements trkf::Surface.

Definition at line 34 of file SurfLine.cxx.

35  {
36  // Get slope parameters and error matrix.
37 
38  double phi = vec(2);
39  double eta = vec(3);
40  double epp = err(2, 2); // sigma^2(phi,phi)
41  double ehh = err(3, 3); // sigma^2(eta,eta)
42  double ehp = err(3, 2); // sigma^2(eta,phi)
43 
44  // Calculate error matrix of pointing unit vector in some coordinate system.
45 
46  double sh = 1. / std::cosh(eta); // sech(eta)
47  double sh2 = sh * sh;
48  double sh3 = sh * sh2;
49  double sh4 = sh * sh3;
50 
51  double th = std::tanh(eta);
52  double th2 = th * th;
53 
54  double cphi = std::cos(phi);
55  double cphi2 = cphi * cphi;
56 
57  double sphi = std::sin(phi);
58  double sphi2 = sphi * sphi;
59 
60  double vxx = sh2 * th2 * cphi2 * ehh + sh2 * sphi2 * epp + 2. * sh2 * th * sphi * cphi * ehp;
61  double vyy = sh2 * th2 * sphi2 * ehh + sh2 * cphi2 * epp - 2. * sh2 * th * sphi * cphi * ehp;
62  double vzz = sh4 * epp;
63 
64  double vxy =
65  sh2 * th2 * sphi * cphi * ehh - sh2 * sphi * cphi * epp - sh2 * th * (cphi2 - sphi2) * ehp;
66  double vyz = -sh3 * th * sphi * ehh + sh3 * cphi * ehp;
67  double vxz = -sh3 * th * cphi * ehh - sh3 * sphi * ehp;
68 
69  // For debugging. The determinant of the error matrix should be zero.
70 
71  // double det = vxx*vyy*vzz + 2.*vxy*vyz*vxz - vxx*vyz*vyz - vyy*vxz*vxz - vzz*vxy*vxy;
72 
73  // Calculate square root of the largest eigenvalue of error matrix.
74 
75  double ddd2 = vxx * vxx + vyy * vyy + vzz * vzz - 2. * vxx * vyy - 2. * vxx * vzz -
76  2. * vyy * vzz + 4. * vxy * vxy + 4. * vyz * vyz + 4. * vxz * vxz;
77  double ddd = sqrt(ddd2 > 0. ? ddd2 : 0.);
78  double lambda2 = 0.5 * (vxx + vyy + vzz + ddd);
79  double lambda = sqrt(lambda2 > 0. ? lambda2 : 0.);
80 
81  return lambda;
82  }
std::ostream & trkf::SurfYZLine::Print ( std::ostream &  out) const
virtual

Printout.

Implements trkf::Surface.

Definition at line 282 of file SurfYZLine.cxx.

References fPhi, fX0, fY0, and fZ0.

Referenced by phi().

283  {
284  out << "SurfYZLine{ x0=" << fX0 << ", y0=" << fY0 << ", z0=" << fZ0 << ", phi=" << fPhi << "}";
285  return out;
286  }
double fX0
X origin.
Definition: SurfYZLine.h:138
double fZ0
Z origin.
Definition: SurfYZLine.h:140
double fY0
Y origin.
Definition: SurfYZLine.h:139
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
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 84 of file SurfYZLine.cxx.

References fPhi, fX0, fY0, and fZ0.

Referenced by getPosition(), and phi().

85  {
86  double sinphi = std::sin(fPhi);
87  double cosphi = std::cos(fPhi);
88 
89  // x = x0 + u
90  xyz[0] = fX0 + uvw[0];
91 
92  // y = y0 + v*cos(phi) - w*sin(phi)
93  xyz[1] = fY0 + uvw[1] * cosphi - uvw[2] * sinphi;
94 
95  // z = z0 + v*sin(phi) + w*cos(phi)
96  xyz[2] = fZ0 + uvw[1] * sinphi + uvw[2] * cosphi;
97  }
double fX0
X origin.
Definition: SurfYZLine.h:138
double fZ0
Z origin.
Definition: SurfYZLine.h:140
double fY0
Y origin.
Definition: SurfYZLine.h:139
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
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 62 of file SurfYZLine.cxx.

References fPhi, fX0, fY0, and fZ0.

Referenced by distanceTo(), and phi().

63  {
64  double sinphi = std::sin(fPhi);
65  double cosphi = std::cos(fPhi);
66 
67  // u = x-x0
68  uvw[0] = xyz[0] - fX0;
69 
70  // v = (y-y0)*cos(phi) + (z-z0)*sin(phi)
71  uvw[1] = (xyz[1] - fY0) * cosphi + (xyz[2] - fZ0) * sinphi;
72 
73  // w = -(y-y0)*sin(phi) + (z-z0)*cos(phi)
74  uvw[2] = -(xyz[1] - fY0) * sinphi + (xyz[2] - fZ0) * cosphi;
75  }
double fX0
X origin.
Definition: SurfYZLine.h:138
double fZ0
Z origin.
Definition: SurfYZLine.h:140
double fY0
Y origin.
Definition: SurfYZLine.h:139
double fPhi
Rotation angle about x-axis.
Definition: SurfYZLine.h:141
double trkf::SurfYZLine::x0 ( ) const
inline

X origin.

Definition at line 90 of file SurfYZLine.h.

References fX0.

Referenced by isEqual(), and trkf::PropYZLine::short_vec_prop().

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

Y origin.

Definition at line 91 of file SurfYZLine.h.

References fY0.

Referenced by isEqual(), and trkf::PropYZLine::short_vec_prop().

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

Z origin.

Definition at line 92 of file SurfYZLine.h.

References fZ0.

Referenced by isEqual(), and trkf::PropYZLine::short_vec_prop().

Member Data Documentation

double trkf::SurfYZLine::fPhi
private

Rotation angle about x-axis.

Definition at line 141 of file SurfYZLine.h.

Referenced by getMomentum(), isEqual(), isParallel(), phi(), Print(), toGlobal(), and toLocal().

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

Phi tolerance for parallel.

Definition at line 133 of file SurfYZLine.h.

Referenced by isEqual(), and isParallel().

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

Separation tolerance for equal.

Definition at line 134 of file SurfYZLine.h.

Referenced by isEqual().

double trkf::SurfYZLine::fX0
private

X origin.

Definition at line 138 of file SurfYZLine.h.

Referenced by isEqual(), Print(), toGlobal(), toLocal(), and x0().

double trkf::SurfYZLine::fY0
private

Y origin.

Definition at line 139 of file SurfYZLine.h.

Referenced by isEqual(), Print(), toGlobal(), toLocal(), and y0().

double trkf::SurfYZLine::fZ0
private

Z origin.

Definition at line 140 of file SurfYZLine.h.

Referenced by isEqual(), Print(), toGlobal(), toLocal(), and z0().


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