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

#include <SurfXYZPlane.h>

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

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...
 
- Public Member Functions inherited from trkf::SurfPlane
 SurfPlane ()
 Default constructor. More...
 
virtual ~SurfPlane ()
 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 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...
 

Additional Inherited Members

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

Detailed Description

Definition at line 83 of file SurfXYZPlane.h.

Constructor & Destructor Documentation

trkf::SurfXYZPlane::SurfXYZPlane ( )

Default constructor.

Definition at line 25 of file SurfXYZPlane.cxx.

25  :
26  fX0(0.),
27  fY0(0.),
28  fZ0(0.),
29  fPhi(0.),
30  fTheta(0.)
31  {}
double fX0
X origin.
Definition: SurfXYZPlane.h:148
double fY0
Y origin.
Definition: SurfXYZPlane.h:149
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
double fZ0
Z origin.
Definition: SurfXYZPlane.h:150
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 41 of file SurfXYZPlane.cxx.

41  :
42  fX0(x0),
43  fY0(y0),
44  fZ0(z0),
45  fPhi(phi),
46  fTheta(theta)
47  {}
double y0() const
Y origin.
Definition: SurfXYZPlane.h:102
double fX0
X origin.
Definition: SurfXYZPlane.h:148
double phi() const
Rot. angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:104
double x0() const
X origin.
Definition: SurfXYZPlane.h:101
double theta() const
Rot. angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:105
double fY0
Y origin.
Definition: SurfXYZPlane.h:149
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
double z0() const
Z origin.
Definition: SurfXYZPlane.h:103
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
double fZ0
Z origin.
Definition: SurfXYZPlane.h:150
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 56 of file SurfXYZPlane.cxx.

57  :
58  fX0(x0),
59  fY0(y0),
60  fZ0(z0),
61  fPhi(0.),
62  fTheta(0.)
63  {
64  // Calculate angles.
65 
66  double nyz = std::hypot(ny, nz);
67  fTheta = atan2(nx, nyz);
68  fPhi = 0.;
69  if(nyz != 0.)
70  fPhi = atan2(-ny, nz);
71  }
double y0() const
Y origin.
Definition: SurfXYZPlane.h:102
std::enable_if_t< std::is_arithmetic_v< T >, T > hypot(T x, T y)
Definition: hypot.h:60
double fX0
X origin.
Definition: SurfXYZPlane.h:148
double x0() const
X origin.
Definition: SurfXYZPlane.h:101
double fY0
Y origin.
Definition: SurfXYZPlane.h:149
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
double z0() const
Z origin.
Definition: SurfXYZPlane.h:103
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
double fZ0
Z origin.
Definition: SurfXYZPlane.h:150
trkf::SurfXYZPlane::~SurfXYZPlane ( )
virtual

Destructor.

Definition at line 74 of file SurfXYZPlane.cxx.

75  {}

Member Function Documentation

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

Clone method.

Implements trkf::Surface.

Definition at line 78 of file SurfXYZPlane.cxx.

79  {
80  return new SurfXYZPlane(*this);
81  }
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 253 of file SurfXYZPlane.cxx.

254  {
255  // Check if the other surface is parallel to this one.
256 
257  bool parallel = isParallel(surf);
258  if(!parallel)
259  throw cet::exception("SurfXYZPlane") << "Attempt to find distance to non-parallel surface.\n";
260 
261  // Find the origin of the other surface in global coordinates,
262  // then convert to our local coordinates.
263 
264  double otheruvw[3] = {0., 0., 0.};
265  double xyz[3];
266  double myuvw[3];
267  surf.toGlobal(otheruvw, xyz);
268  toLocal(xyz, myuvw);
269 
270  // Distance is local w-coordinate of other surface origin.
271 
272  return myuvw[2];
273  }
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
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 167 of file SurfXYZPlane.cxx.

169  {
170 
171  // Get momentum.
172 
173  double invp = std::abs(vec(4));
174  double p = 1. / std::max(invp, 1.e-3); // Capped at 1000. GeV/c.
175 
176  // Get track slope parameters.
177 
178  double dudw = vec(2);
179  double dvdw = vec(3);
180 
181  // Calculate dw/ds.
182 
183  double dwds = 1. / std::sqrt(1. + dudw*dudw + dvdw*dvdw);
184  TrackDirection realdir = getDirection(vec, dir); // Should be same as original direction.
185  if(realdir == BACKWARD)
186  dwds = -dwds;
187  else if(realdir != FORWARD)
188  throw cet::exception("SurfXYZPlane") << "Track direction not specified.\n";
189 
190  // Calculate momentum vector in local coordinate system.
191 
192  double pu = p * dudw * dwds;
193  double pv = p * dvdw * dwds;
194  double pw = p * dwds;
195 
196  // Rotate momentum to global coordinte system.
197 
198  double sinth = std::sin(fTheta);
199  double costh = std::cos(fTheta);
200  double sinphi = std::sin(fPhi);
201  double cosphi = std::cos(fPhi);
202 
203  mom[0] = pu*costh + pw*sinth;
204  mom[1] = pu*sinth*sinphi + pv*cosphi - pw*costh*sinphi;
205  mom[2] = -pu*sinth*cosphi + pv*sinphi + pw*costh*cosphi;
206 
207  return;
208  }
TrackDirection
Track direction enum.
Definition: Surface.h:56
string dir
T abs(T value)
const double e
p
Definition: test.py:223
static int max(int a, int b)
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
virtual TrackDirection getDirection(const TrackVector &, TrackDirection dir=UNKNOWN) const
Get direction of track (default UNKNOWN).
Definition: Surface.h:85
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
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 144 of file SurfXYZPlane.cxx.

145  {
146  // Get position in local coordinate system.
147 
148  double uvw[3];
149  uvw[0] = vec(0);
150  uvw[1] = vec(1);
151  uvw[2] = 0.;
152 
153  // Transform to global coordinate system.
154 
155  toGlobal(uvw, xyz);
156  return;
157  }
virtual void toGlobal(const double uvw[3], double xyz[3]) const
Transform local to global coordinates.
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 286 of file SurfXYZPlane.cxx.

287  {
288  bool result = false;
289 
290  // Test if the other surface is parallel.
291 
292  bool parallel = isParallel(surf);
293  if(parallel) {
294  double dist = distanceTo(surf);
295  if(std::abs(dist) <= fSepTolerance)
296  result = true;
297  }
298 
299  return result;
300  }
static QCString result
virtual bool isParallel(const Surface &surf) const
Test whether two surfaces are parallel, within tolerance.
T abs(T value)
static double fSepTolerance
Separation tolerance for equal.
Definition: SurfXYZPlane.h:144
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 220 of file SurfXYZPlane.cxx.

221  {
222  bool result = false;
223 
224  // Test if the other surface is a SurfXYZPlane.
225 
226  const SurfXYZPlane* psurf = dynamic_cast<const SurfXYZPlane*>(&surf);
227  if(psurf != 0) {
228 
229  // Test whether surface angle parameters are the same
230  // with tolerance.
231 
232  double delta_phi = TVector2::Phi_mpi_pi(fPhi - psurf->phi());
233  double delta_theta = fTheta - psurf->theta();
234  if(std::abs(delta_phi) <= fPhiTolerance && std::abs(delta_theta) <= fThetaTolerance)
235  result = true;
236  }
237  return result;
238  }
static QCString result
static double fPhiTolerance
Phi tolerance for parallel.
Definition: SurfXYZPlane.h:142
T abs(T value)
SurfXYZPlane()
Default constructor.
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
static double fThetaTolerance
Theta tolerance for parallel.
Definition: SurfXYZPlane.h:143
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
bool trkf::SurfXYZPlane::isTrackValid ( const TrackVector vec) const
virtual

Surface-specific tests of validity of track parameters.

Implements trkf::Surface.

Definition at line 84 of file SurfXYZPlane.cxx.

85  {
86  return true;
87  }
double trkf::SurfXYZPlane::phi ( ) const
inline

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

Definition at line 104 of file SurfXYZPlane.h.

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

Printout.

Implements trkf::Surface.

Definition at line 303 of file SurfXYZPlane.cxx.

304  {
305  out << "SurfXYZPlane{ x0=" << fX0 << ", y0=" << fY0 << ", z0=" << fZ0
306  << ", phi=" << fPhi << ", theta=" << fTheta << "}";
307  return out;
308  }
double fX0
X origin.
Definition: SurfXYZPlane.h:148
double fY0
Y origin.
Definition: SurfXYZPlane.h:149
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
double fZ0
Z origin.
Definition: SurfXYZPlane.h:150
double trkf::SurfXYZPlane::theta ( ) const
inline

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

Definition at line 105 of file SurfXYZPlane.h.

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 120 of file SurfXYZPlane.cxx.

121  {
122  double sinth = std::sin(fTheta);
123  double costh = std::cos(fTheta);
124  double sinphi = std::sin(fPhi);
125  double cosphi = std::cos(fPhi);
126 
127  // x = x0 + u*cos(theta) + w*sin(theta)
128  xyz[0] = fX0 + uvw[0]*costh + uvw[2]*sinth;
129 
130  // y = y0 + u*sin(theta)*sin(phi) + v*cos(phi) - w*cos(theta)*sin(phi)
131  xyz[1] = fY0 + uvw[0]*sinth*sinphi + uvw[1]*cosphi - uvw[2]*costh*sinphi;
132 
133  // z = z0 - u*sin(theta)*cos(phi) + v*sin(phi) + w*cos(theta)*cos(phi)
134  xyz[2] = fZ0 - uvw[0]*sinth*cosphi + uvw[1]*sinphi + uvw[2]*costh*cosphi;
135  }
double fX0
X origin.
Definition: SurfXYZPlane.h:148
double fY0
Y origin.
Definition: SurfXYZPlane.h:149
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
double fZ0
Z origin.
Definition: SurfXYZPlane.h:150
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 96 of file SurfXYZPlane.cxx.

97  {
98  double sinth = std::sin(fTheta);
99  double costh = std::cos(fTheta);
100  double sinphi = std::sin(fPhi);
101  double cosphi = std::cos(fPhi);
102 
103  // u = (x-x0)*cos(theta) + (y-y0)*sin(theta)*sin(phi) - (z-z0)*sin(theta)*cos(phi)
104  uvw[0] = (xyz[0]-fX0)*costh + (xyz[1]-fY0)*sinth*sinphi - (xyz[2]-fZ0)*sinth*cosphi;
105 
106  // v = (y-y0)*cos(phi) + (z-z0)*sin(phi)
107  uvw[1] = (xyz[1]-fY0)*cosphi + (xyz[2]-fZ0)*sinphi;
108 
109  // w = (x-x0)*sin(theta) - (y-y0)*cos(theta)*sin(phi) + (z-z0)*cos(theta)*cos(phi)
110  uvw[2] = (xyz[0]-fX0)*sinth - (xyz[1]-fY0)*costh*sinphi + (xyz[2]-fZ0)*costh*cosphi;
111  }
double fX0
X origin.
Definition: SurfXYZPlane.h:148
double fY0
Y origin.
Definition: SurfXYZPlane.h:149
double fTheta
Rotation angle about y&#39;-axis (projected Lorentz angle).
Definition: SurfXYZPlane.h:152
double fPhi
Rotation angle about x-axis (wire angle).
Definition: SurfXYZPlane.h:151
double fZ0
Z origin.
Definition: SurfXYZPlane.h:150
double trkf::SurfXYZPlane::x0 ( ) const
inline

X origin.

Definition at line 101 of file SurfXYZPlane.h.

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

Y origin.

Definition at line 102 of file SurfXYZPlane.h.

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

Z origin.

Definition at line 103 of file SurfXYZPlane.h.

Member Data Documentation

double trkf::SurfXYZPlane::fPhi
private

Rotation angle about x-axis (wire angle).

Definition at line 151 of file SurfXYZPlane.h.

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

Phi tolerance for parallel.

Definition at line 142 of file SurfXYZPlane.h.

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

Separation tolerance for equal.

Definition at line 144 of file SurfXYZPlane.h.

double trkf::SurfXYZPlane::fTheta
private

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

Definition at line 152 of file SurfXYZPlane.h.

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

Theta tolerance for parallel.

Definition at line 143 of file SurfXYZPlane.h.

double trkf::SurfXYZPlane::fX0
private

X origin.

Definition at line 148 of file SurfXYZPlane.h.

double trkf::SurfXYZPlane::fY0
private

Y origin.

Definition at line 149 of file SurfXYZPlane.h.

double trkf::SurfXYZPlane::fZ0
private

Z origin.

Definition at line 150 of file SurfXYZPlane.h.


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