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

#include "SurfXYZPlane.h"

Inheritance diagram for trkf::SurfXYZPlane:
trkf::SurfPlane trkf::Surface

Public Types

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

Public Member Functions

 SurfXYZPlane ()
 Default constructor. More...
 
 SurfXYZPlane (double x0, double y0, double z0, double phi, double theta)
 Initializing constructor (angles). More...
 
 SurfXYZPlane (double x0, double y0, double z0, double nx, double ny, double nz)
 Initializing constructor (normal vector). More...
 
virtual ~SurfXYZPlane ()
 Destructor. More...
 
double x0 () const
 X origin. More...
 
double y0 () const
 Y origin. More...
 
double z0 () const
 Z origin. More...
 
double phi () const
 Rot. angle about x-axis (wire angle). More...
 
double theta () const
 Rot. angle about y'-axis (projected Lorentz angle). 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 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 forward 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 TrackVector getDiff (const TrackVector &vec1, const TrackVector &vec2) const
 Calculate difference of two track parameter vectors. 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 (wire angle). More...
 
double fTheta
 Rotation angle about y'-axis (projected Lorentz angle). More...
 

Static Private Attributes

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

Detailed Description

Definition at line 83 of file SurfXYZPlane.h.

Member Enumeration Documentation

Track direction enum.

Enumerator
FORWARD 
BACKWARD 
UNKNOWN 

Definition at line 54 of file Surface.h.

Constructor & Destructor Documentation

trkf::SurfXYZPlane::SurfXYZPlane ( )

Default constructor.

Definition at line 25 of file SurfXYZPlane.cxx.

Referenced by clone().

25 : fX0(0.), fY0(0.), fZ0(0.), fPhi(0.), fTheta(0.) {}
double fX0
X origin.
Definition: SurfXYZPlane.h:145
double fY0
Y origin.
Definition: SurfXYZPlane.h:146
double fTheta
Rotation angle about y'-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
double fZ0
Z origin.
Definition: SurfXYZPlane.h:147
trkf::SurfXYZPlane::SurfXYZPlane ( double  x0,
double  y0,
double  z0,
double  phi,
double  theta 
)

Initializing constructor (angles).

Initializing constructor.

Arguments:

x0, y0, z0 - Global coordinates of local origin. phi - Rotation angle about x-axis (wire angle). theta - Rotation angle about y'-axis (projected Lorentz angle).

Definition at line 35 of file SurfXYZPlane.cxx.

36  : fX0(x0), fY0(y0), fZ0(z0), fPhi(phi), fTheta(theta)
37  {}
double y0() const
Y origin.
Definition: SurfXYZPlane.h:99
double fX0
X origin.
Definition: SurfXYZPlane.h:145
double phi() const
Rot. angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:101
double x0() const
X origin.
Definition: SurfXYZPlane.h:98
double theta() const
Rot. angle about y'-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:102
double fY0
Y origin.
Definition: SurfXYZPlane.h:146
double fTheta
Rotation angle about y'-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
double z0() const
Z origin.
Definition: SurfXYZPlane.h:100
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
double fZ0
Z origin.
Definition: SurfXYZPlane.h:147
trkf::SurfXYZPlane::SurfXYZPlane ( double  x0,
double  y0,
double  z0,
double  nx,
double  ny,
double  nz 
)

Initializing constructor (normal vector).

Initializing constructor (normal vector).

Arguments:

x0, y0, z0 - Global coordinates of local origin. nx, ny, nz - Normal vector in global coordinate system.

Definition at line 46 of file SurfXYZPlane.cxx.

References fPhi, and fTheta.

47  : fX0(x0), fY0(y0), fZ0(z0), fPhi(0.), fTheta(0.)
48  {
49  // Calculate angles.
50 
51  double nyz = std::hypot(ny, nz);
52  fTheta = atan2(nx, nyz);
53  fPhi = 0.;
54  if (nyz != 0.) fPhi = atan2(-ny, nz);
55  }
double y0() const
Y origin.
Definition: SurfXYZPlane.h:99
double fX0
X origin.
Definition: SurfXYZPlane.h:145
double x0() const
X origin.
Definition: SurfXYZPlane.h:98
double fY0
Y origin.
Definition: SurfXYZPlane.h:146
double fTheta
Rotation angle about y'-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
double z0() const
Z origin.
Definition: SurfXYZPlane.h:100
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
double fZ0
Z origin.
Definition: SurfXYZPlane.h:147
trkf::SurfXYZPlane::~SurfXYZPlane ( )
virtual

Destructor.

Definition at line 58 of file SurfXYZPlane.cxx.

58 {}

Member Function Documentation

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

Clone method.

Implements trkf::Surface.

Definition at line 61 of file SurfXYZPlane.cxx.

References SurfXYZPlane().

Referenced by theta().

62  {
63  return new SurfXYZPlane(*this);
64  }
SurfXYZPlane()
Default constructor.
double trkf::SurfXYZPlane::distanceTo ( const Surface surf) const
virtual

Find perpendicular forward distance to a parallel surface.

Find perpendicular forward distance to a parallel surface.

Throw an exception if the other surface is not parallel.

Assuming the other surface is parallel, the distance is simply the w-coordinate of the other surface, and is signed.

Arguments:

surf - Other surface.

Returned value: Distance.

Implements trkf::Surface.

Definition at line 237 of file SurfXYZPlane.cxx.

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

Referenced by isEqual(), and theta().

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("SurfXYZPlane") << "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 is local w-coordinate of other surface origin.
255 
256  return myuvw[2];
257  }
virtual bool isParallel(const Surface &surf) const
Test whether two surfaces are parallel, within tolerance.
virtual void toLocal(const double xyz[3], double uvw[3]) const
Transform global to local coordinates.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
TrackVector trkf::Surface::getDiff ( const TrackVector vec1,
const TrackVector vec2 
) const
virtualinherited

Calculate difference of two track parameter vectors.

Calculate difference of two track parameter vectors. This method has a default implementation which is just the numeric difference. Surfaces that require a more sophisticated difference (e.g. phi-wrap difference) should override this method.

Reimplemented in trkf::SurfYZLine.

Definition at line 31 of file Surface.cxx.

32  {
33  return vec1 - vec2;
34  }
virtual TrackDirection trkf::Surface::getDirection ( const TrackVector ,
TrackDirection  dir = UNKNOWN 
) const
inlinevirtualinherited
void trkf::SurfXYZPlane::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.

Implements trkf::Surface.

Definition at line 152 of file SurfXYZPlane.cxx.

References util::abs(), trkf::Surface::BACKWARD, e, trkf::Surface::FORWARD, fPhi, fTheta, and trkf::Surface::getDirection().

Referenced by theta().

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 slope parameters.
161 
162  double dudw = vec(2);
163  double dvdw = vec(3);
164 
165  // Calculate dw/ds.
166 
167  double dwds = 1. / std::sqrt(1. + dudw * dudw + dvdw * dvdw);
168  TrackDirection realdir = getDirection(vec, dir); // Should be same as original direction.
169  if (realdir == BACKWARD)
170  dwds = -dwds;
171  else if (realdir != FORWARD)
172  throw cet::exception("SurfXYZPlane") << "Track direction not specified.\n";
173 
174  // Calculate momentum vector in local coordinate system.
175 
176  double pu = p * dudw * dwds;
177  double pv = p * dvdw * dwds;
178  double pw = p * dwds;
179 
180  // Rotate momentum to global coordinte system.
181 
182  double sinth = std::sin(fTheta);
183  double costh = std::cos(fTheta);
184  double sinphi = std::sin(fPhi);
185  double cosphi = std::cos(fPhi);
186 
187  mom[0] = pu * costh + pw * sinth;
188  mom[1] = pu * sinth * sinphi + pv * cosphi - pw * costh * sinphi;
189  mom[2] = -pu * sinth * cosphi + pv * sinphi + pw * costh * cosphi;
190 
191  return;
192  }
TrackDirection
Track direction enum.
Definition: Surface.h:54
constexpr auto abs(T v)
Returns the absolute value of the argument.
TDirectory * dir
Definition: macro.C:5
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
virtual TrackDirection getDirection(const TrackVector &, TrackDirection dir=UNKNOWN) const
Get direction of track (default UNKNOWN).
Definition: Surface.h:83
Float_t e
Definition: plot.C:35
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
void trkf::SurfXYZPlane::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 129 of file SurfXYZPlane.cxx.

References toGlobal().

Referenced by theta().

130  {
131  // Get position in local coordinate system.
132 
133  double uvw[3];
134  uvw[0] = vec(0);
135  uvw[1] = vec(1);
136  uvw[2] = 0.;
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.
void trkf::SurfPlane::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 83 of file SurfPlane.cxx.

84  {
85  err.resize(5, false);
86  err.clear();
87  err(0, 0) = 1000.;
88  err(1, 1) = 1000.;
89  err(2, 2) = 0.25;
90  err(3, 3) = 0.25;
91  err(4, 4) = 10.;
92  }
bool trkf::SurfXYZPlane::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 parallel and having zero separation, within tolerance. Note that this definition of equality allows the two surfaces to have different origins.

Arguments:

surf - Other surface.

Returned values: true if equal.

Implements trkf::Surface.

Definition at line 270 of file SurfXYZPlane.cxx.

References util::abs(), larg4::dist(), distanceTo(), fSepTolerance, and isParallel().

Referenced by theta().

271  {
272  bool result = false;
273 
274  // Test if the other surface is parallel.
275 
276  bool parallel = isParallel(surf);
277  if (parallel) {
278  double dist = distanceTo(surf);
279  if (std::abs(dist) <= fSepTolerance) result = true;
280  }
281 
282  return result;
283  }
constexpr auto abs(T v)
Returns the absolute value of the argument.
virtual bool isParallel(const Surface &surf) const
Test whether two surfaces are parallel, within tolerance.
static double fSepTolerance
Separation tolerance for equal.
Definition: SurfXYZPlane.h:141
virtual double distanceTo(const Surface &surf) const
Find perpendicular forward distance to a parallel surface.
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
bool trkf::SurfXYZPlane::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 SurfXYZPlane.

Arguments:

surf - Other surface.

Returned value: true if parallel.

Implements trkf::Surface.

Definition at line 204 of file SurfXYZPlane.cxx.

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

Referenced by distanceTo(), isEqual(), and theta().

205  {
206  bool result = false;
207 
208  // Test if the other surface is a SurfXYZPlane.
209 
210  const SurfXYZPlane* psurf = dynamic_cast<const SurfXYZPlane*>(&surf);
211  if (psurf != 0) {
212 
213  // Test whether surface angle parameters are the same
214  // with tolerance.
215 
216  double delta_phi = TVector2::Phi_mpi_pi(fPhi - psurf->phi());
217  double delta_theta = fTheta - psurf->theta();
218  if (std::abs(delta_phi) <= fPhiTolerance && std::abs(delta_theta) <= fThetaTolerance)
219  result = true;
220  }
221  return result;
222  }
constexpr auto abs(T v)
Returns the absolute value of the argument.
static double fPhiTolerance
Phi tolerance for parallel.
Definition: SurfXYZPlane.h:139
SurfXYZPlane()
Default constructor.
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
static double fThetaTolerance
Theta tolerance for parallel.
Definition: SurfXYZPlane.h:140
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
bool trkf::SurfXYZPlane::isTrackValid ( const TrackVector vec) const
virtual

Surface-specific tests of validity of track parameters.

Implements trkf::Surface.

Definition at line 67 of file SurfXYZPlane.cxx.

Referenced by theta().

68  {
69  return true;
70  }
double trkf::SurfXYZPlane::phi ( ) const
inline

Rot. angle about x-axis (wire angle).

Definition at line 101 of file SurfXYZPlane.h.

References fPhi.

Referenced by isParallel(), trkf::PropXYZPlane::origin_vec_prop(), and trkf::PropXYZPlane::short_vec_prop().

double trkf::SurfPlane::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 SurfPlane.cxx.

References den.

35  {
36  // Get slope parameters and error matrix.
37 
38  double xp = vec(2);
39  double yp = vec(3);
40  double exx = err(2, 2);
41  double eyy = err(3, 3);
42  double exy = err(3, 2);
43 
44  // Calculate error matrix of pointing unit vector in some coordinate system.
45 
46  double den = 1. + xp * xp + yp * yp;
47  double den3 = den * den * den;
48 
49  double vxx = ((1. + yp * yp) * (1. + yp * yp) * exx + xp * xp * yp * yp * eyy -
50  2. * xp * yp * (1. + yp * yp) * exy) /
51  den3;
52  double vyy = (xp * xp * yp * yp * exx + (1. + xp * xp) * (1. + xp * xp) * eyy -
53  2. * xp * yp * (1. + xp * xp) * exy) /
54  den3;
55  double vzz = (xp * xp * exx + yp * yp * eyy + 2. * xp * yp * exy) / den3;
56 
57  double vxy = (-xp * yp * (1. + yp * yp) * exx - xp * yp * (1. + xp * xp) * eyy +
58  (1. + xp * xp + yp * yp + 2. * xp * xp * yp * yp) * exy) /
59  den3;
60  double vyz =
61  (xp * xp * yp * exx - yp * (1. + xp * xp) * eyy - xp * (1. + xp * xp - yp * yp) * exy) / den3;
62  double vxz =
63  (-xp * (1. + yp * yp) * exx + xp * yp * yp * eyy - yp * (1. - xp * xp + yp * yp) * exy) /
64  den3;
65 
66  // Calculate square root of the largest eigenvalue of error matrix.
67 
68  double ddd2 = vxx * vxx + vyy * vyy + vzz * vzz - 2. * vxx * vyy - 2. * vxx * vzz -
69  2. * vyy * vzz + 4. * vxy * vxy + 4. * vyz * vyz + 4. * vxz * vxz;
70  double ddd = sqrt(ddd2 > 0. ? ddd2 : 0.);
71  double lambda2 = 0.5 * (vxx + vyy + vzz + ddd);
72  double lambda = sqrt(lambda2 > 0. ? lambda2 : 0.);
73 
74  return lambda;
75  }
Float_t den
Definition: plot.C:35
std::ostream & trkf::SurfXYZPlane::Print ( std::ostream &  out) const
virtual

Printout.

Implements trkf::Surface.

Definition at line 286 of file SurfXYZPlane.cxx.

References fPhi, fTheta, fX0, fY0, and fZ0.

Referenced by theta().

287  {
288  out << "SurfXYZPlane{ x0=" << fX0 << ", y0=" << fY0 << ", z0=" << fZ0 << ", phi=" << fPhi
289  << ", theta=" << fTheta << "}";
290  return out;
291  }
double fX0
X origin.
Definition: SurfXYZPlane.h:145
double fY0
Y origin.
Definition: SurfXYZPlane.h:146
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
double fZ0
Z origin.
Definition: SurfXYZPlane.h:147
double trkf::SurfXYZPlane::theta ( ) const
inline
void trkf::SurfXYZPlane::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 105 of file SurfXYZPlane.cxx.

References fPhi, fTheta, fX0, fY0, and fZ0.

Referenced by getPosition(), and theta().

106  {
107  double sinth = std::sin(fTheta);
108  double costh = std::cos(fTheta);
109  double sinphi = std::sin(fPhi);
110  double cosphi = std::cos(fPhi);
111 
112  // x = x0 + u*cos(theta) + w*sin(theta)
113  xyz[0] = fX0 + uvw[0] * costh + uvw[2] * sinth;
114 
115  // y = y0 + u*sin(theta)*sin(phi) + v*cos(phi) - w*cos(theta)*sin(phi)
116  xyz[1] = fY0 + uvw[0] * sinth * sinphi + uvw[1] * cosphi - uvw[2] * costh * sinphi;
117 
118  // z = z0 - u*sin(theta)*cos(phi) + v*sin(phi) + w*cos(theta)*cos(phi)
119  xyz[2] = fZ0 - uvw[0] * sinth * cosphi + uvw[1] * sinphi + uvw[2] * costh * cosphi;
120  }
double fX0
X origin.
Definition: SurfXYZPlane.h:145
double fY0
Y origin.
Definition: SurfXYZPlane.h:146
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
double fZ0
Z origin.
Definition: SurfXYZPlane.h:147
void trkf::SurfXYZPlane::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 79 of file SurfXYZPlane.cxx.

References fPhi, fTheta, fX0, fY0, and fZ0.

Referenced by distanceTo(), and theta().

80  {
81  double sinth = std::sin(fTheta);
82  double costh = std::cos(fTheta);
83  double sinphi = std::sin(fPhi);
84  double cosphi = std::cos(fPhi);
85 
86  // u = (x-x0)*cos(theta) + (y-y0)*sin(theta)*sin(phi) - (z-z0)*sin(theta)*cos(phi)
87  uvw[0] =
88  (xyz[0] - fX0) * costh + (xyz[1] - fY0) * sinth * sinphi - (xyz[2] - fZ0) * sinth * cosphi;
89 
90  // v = (y-y0)*cos(phi) + (z-z0)*sin(phi)
91  uvw[1] = (xyz[1] - fY0) * cosphi + (xyz[2] - fZ0) * sinphi;
92 
93  // w = (x-x0)*sin(theta) - (y-y0)*cos(theta)*sin(phi) + (z-z0)*cos(theta)*cos(phi)
94  uvw[2] =
95  (xyz[0] - fX0) * sinth - (xyz[1] - fY0) * costh * sinphi + (xyz[2] - fZ0) * costh * cosphi;
96  }
double fX0
X origin.
Definition: SurfXYZPlane.h:145
double fY0
Y origin.
Definition: SurfXYZPlane.h:146
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:149
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:148
double fZ0
Z origin.
Definition: SurfXYZPlane.h:147
double trkf::SurfXYZPlane::x0 ( ) const
inline

X origin.

Definition at line 98 of file SurfXYZPlane.h.

References fX0.

Referenced by trkf::PropXYZPlane::short_vec_prop().

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

Y origin.

Definition at line 99 of file SurfXYZPlane.h.

References fY0.

Referenced by trkf::PropXYZPlane::short_vec_prop().

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

Z origin.

Definition at line 100 of file SurfXYZPlane.h.

References fZ0.

Referenced by trkf::PropXYZPlane::short_vec_prop().

Member Data Documentation

double trkf::SurfXYZPlane::fPhi
private

Rotation angle about x-axis (wire angle).

Definition at line 148 of file SurfXYZPlane.h.

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

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

Phi tolerance for parallel.

Definition at line 139 of file SurfXYZPlane.h.

Referenced by isParallel().

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

Separation tolerance for equal.

Definition at line 141 of file SurfXYZPlane.h.

Referenced by isEqual().

double trkf::SurfXYZPlane::fTheta
private

Rotation angle about y'-axis (projected Lorentz angle).

Definition at line 149 of file SurfXYZPlane.h.

Referenced by getMomentum(), isParallel(), Print(), SurfXYZPlane(), theta(), toGlobal(), and toLocal().

double trkf::SurfXYZPlane::fThetaTolerance = 1.e-10
staticprivate

Theta tolerance for parallel.

Definition at line 140 of file SurfXYZPlane.h.

Referenced by isParallel().

double trkf::SurfXYZPlane::fX0
private

X origin.

Definition at line 145 of file SurfXYZPlane.h.

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

double trkf::SurfXYZPlane::fY0
private

Y origin.

Definition at line 146 of file SurfXYZPlane.h.

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

double trkf::SurfXYZPlane::fZ0
private

Z origin.

Definition at line 147 of file SurfXYZPlane.h.

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


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