All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Trajectory_test.cc
Go to the documentation of this file.
1 /**
2  * @file Trajectory_test.cc
3  * @brief Simple test on a recob::Trajectory object
4  * @author Gianluca Petrillo (petrillo@fnal.gov)
5  * @date January 3, 2016
6  * @version 1.0
7  *
8  * This test simply creates recob::Trajectory objects and verifies that the
9  * values it can access are the right ones.
10  *
11  * See http://www.boost.org/libs/test for the Boost test library home page.
12  *
13  */
14 
15 
16 // Boost libraries
17 /*
18  * Boost: define the name of the module;
19  * and do that before the inclusion of Boost unit test headers
20  * because it will change what they provide.
21  * Among the those, there is a main() function and some wrapping catching
22  * unhandled exceptions and considering them test failures, and probably more.
23  * This also makes fairly complicate to receive parameters from the command line
24  * (for example, a random seed).
25  */
26 #define BOOST_TEST_MODULE ( tracktrajectory_test )
27 #include "boost/test/unit_test.hpp"
28 
29 
30 // LArSoft libraries
33 
34 // ROOT libraries (for the legacy interface)
35 #include "TVector3.h"
36 #include "TMatrixD.h"
37 
38 // C/C++ standard libraries
39 #include <algorithm> // std::equal(), std::accumulate()
40 #include <tuple> // std::tie()
41 #include <iostream>
42 
43 
44 //------------------------------------------------------------------------------
45 //--- Test code
46 //
47 struct Expected_t {
48 
51  bool hasMomenta;
52  double length;
53  double theta;
54  double phi;
55  double zenith;
56  double azimuth;
57 
58 }; // struct Expected_t
59 
60 template <typename T>
61 void CheckValue(T v, T exp, T tol, std::string tag = "") {
62  if (!tag.empty()) BOOST_TEST_MESSAGE(tag);
63  if (std::abs(exp) < (tol / 100.)) BOOST_CHECK_SMALL(v, tol);
64  else BOOST_CHECK_CLOSE(v, exp, tol);
65 } // CheckValue()
66 
67 template <typename VectA, typename VectB>
68 void CheckVectorsEqual(VectA const& v, VectB const& exp) {
69  BOOST_TEST(v.X() == exp.X());
70  BOOST_TEST(v.Y() == exp.Y());
71  BOOST_TEST(v.Z() == exp.Z());
72 } // CheckVectorsEqual()
73 
74 template <typename VectA, typename VectB>
75 void CheckVectorsClose(VectA const& v, VectB const& exp, double tol = 0.01) {
76  CheckValue(v.X(), exp.X(), tol, " X()");
77  CheckValue(v.Y(), exp.Y(), tol, " Y()");
78  CheckValue(v.Z(), exp.Z(), tol, " Z()");
79 } // CheckVectorsClose()
80 
81 
84  m.SetRotationMatrix(Tm);
85  return m;
86 } // makeRotationMatrix()
87 
88 
90  recob::Trajectory const& traj,
91  Expected_t const& expected
92  )
93 {
94 
95  //----------------------------------------------------------------------------
96  const size_t NPoints = traj.NumberTrajectoryPoints();
97  BOOST_TEST(NPoints == expected.positions.size());
98  BOOST_TEST(traj.NumberTrajectoryPoints() == expected.positions.size());
99 
100  for (size_t i = 0; i <= NPoints + 1; ++i) {
101  BOOST_TEST_MESSAGE("HasPoint() position #" << i);
102  BOOST_TEST(traj.HasPoint(i) == (i < NPoints));
103  } // for
104 
105  if (NPoints == 0) return; // nothing else is defined
106 
107  BOOST_TEST(traj.HasMomentum() == expected.hasMomenta);
108 
109  //----------------------------------------------------------------------------
110  BOOST_TEST(traj.LastPoint() == NPoints - 1);
111 
112  //----------------------------------------------------------------------------
113  BOOST_TEST_MESSAGE("Vertex()");
114  CheckVectorsEqual(traj.Vertex(), expected.positions[0]);
115 
116  BOOST_TEST_MESSAGE("Start()");
117  CheckVectorsEqual(traj.Start(), expected.positions[0]);
118 
119  BOOST_TEST_MESSAGE("End()");
120  CheckVectorsEqual(traj.End(), expected.positions[NPoints - 1]);
121 
122  //----------------------------------------------------------------------------
123  for (size_t i = 0; i < NPoints; ++i) {
124 
125  BOOST_TEST_MESSAGE("LocationAtPoint() position #" << i);
126  CheckVectorsEqual(traj.LocationAtPoint(i), expected.positions[i]);
127 
128  } // for
129 
130 
131  //----------------------------------------------------------------------------
132  TVector3 Vstart, Vend;
133  std::tie(Vstart, Vend) = traj.Extent<TVector3>();
134  BOOST_TEST(Vstart[0] == expected.positions[0].X());
135  BOOST_TEST(Vstart[1] == expected.positions[0].Y());
136  BOOST_TEST(Vstart[2] == expected.positions[0].Z());
137  BOOST_TEST(Vend[0] == expected.positions[NPoints - 1].X());
138  BOOST_TEST(Vend[1] == expected.positions[NPoints - 1].Y());
139  BOOST_TEST(Vend[2] == expected.positions[NPoints - 1].Z());
140 
141 
143  std::tie(start, end) = traj.Extent(); // assign both start and end
144  BOOST_TEST_MESSAGE("Extent() start");
145  CheckVectorsEqual(start, expected.positions[0]);
146  BOOST_TEST_MESSAGE("Extent() end");
147  CheckVectorsEqual(end, expected.positions[NPoints - 1]);
148 
149 
150  //----------------------------------------------------------------------------
151  BOOST_CHECK_CLOSE(traj.Length(), expected.length, 0.01); // 0.01%
152  if (NPoints >= 2){
153  BOOST_CHECK_CLOSE(traj.Length(1),
154  expected.length - (expected.positions[1] - expected.positions[0]).R(),
155  0.01
156  );
157  } // if
158 
159 
160  //----------------------------------------------------------------------------
161  BOOST_TEST_MESSAGE("VertexDirection()");
162  CheckVectorsClose(traj.VertexDirection(), expected.momenta[0].Unit());
163  BOOST_CHECK_CLOSE(traj.VertexDirection().Mag2(), 1.0, 0.01);
164 
165  BOOST_TEST_MESSAGE("StartDirection()");
166  CheckVectorsClose(traj.StartDirection(), expected.momenta[0].Unit());
167  BOOST_CHECK_CLOSE(traj.StartDirection().Mag2(), 1.0, 0.01);
168 
169  BOOST_TEST_MESSAGE("EndDirection()");
170  CheckVectorsClose(traj.EndDirection(), expected.momenta[NPoints - 1].Unit());
171  BOOST_CHECK_CLOSE(traj.EndDirection().Mag2(), 1.0, 0.01);
172 
173 
174  //----------------------------------------------------------------------------
175  BOOST_CHECK_CLOSE(traj.Theta(), expected.theta, 0.01);
176  BOOST_CHECK_CLOSE(traj.Phi(), expected.phi, 0.01);
177  BOOST_CHECK_CLOSE(traj.ZenithAngle(), expected.zenith, 0.01);
178  BOOST_CHECK_CLOSE(traj.AzimuthAngle(), expected.azimuth, 0.01);
179 
180 
181  //----------------------------------------------------------------------------
182 
183  BOOST_TEST_MESSAGE("VertexMomentumVector()");
184  CheckVectorsClose(traj.VertexMomentumVector(), expected.momenta[0]);
185 
186  BOOST_TEST_MESSAGE("StartMomentumVector()");
187  CheckVectorsClose(traj.StartMomentumVector(), expected.momenta[0]);
188 
189  BOOST_TEST_MESSAGE("EndMomentumVector()");
190  CheckVectorsClose(traj.EndMomentumVector(), expected.momenta[NPoints - 1]);
191 
192 
193  //----------------------------------------------------------------------------
194  BOOST_TEST_MESSAGE("VertexMomentum()");
195  BOOST_CHECK_CLOSE(traj.VertexMomentum(), expected.momenta[0].R(), 0.01);
196 
197  BOOST_TEST_MESSAGE("StartMomentum()");
198  BOOST_CHECK_CLOSE(traj.StartMomentum(), expected.momenta[0].R(), 0.01);
199 
200  BOOST_TEST_MESSAGE("EndMomentum()");
201  BOOST_CHECK_CLOSE
202  (traj.EndMomentum(), expected.momenta[NPoints - 1].R(), 0.01);
203 
204 
205  //----------------------------------------------------------------------------
206  for (size_t i = 0; i < NPoints; ++i) {
207 
208  BOOST_TEST_MESSAGE("DirectionAtPoint() position #" << i);
209  CheckVectorsClose(traj.DirectionAtPoint(i), expected.momenta[i].Unit());
210 
211  } // for
212 
213 
214  //----------------------------------------------------------------------------
215  for (size_t i = 0; i < NPoints; ++i) {
216 
217  BOOST_TEST_MESSAGE("MomentumVectorAtPoint() position #" << i);
218  if (traj.HasMomentum())
219  CheckVectorsClose(traj.MomentumVectorAtPoint(i), expected.momenta[i]);
220  else {
222  (traj.MomentumVectorAtPoint(i), expected.momenta[i].Unit());
223  }
224 
225  BOOST_TEST_MESSAGE("MomentumAtPoint() position #" << i);
226  if (traj.HasMomentum()) {
227  BOOST_CHECK_CLOSE
228  (traj.MomentumAtPoint(i), expected.momenta[i].R(), 0.01);
229  }
230  else
231  BOOST_CHECK_CLOSE(traj.MomentumAtPoint(i), 1.0, 0.01);
232 
233  } // for
234 
235 
236  //----------------------------------------------------------------------------
237  TVector3 AstartDir, AendDir;
238  std::tie(AstartDir, AendDir) = traj.Direction<TVector3>();
239  BOOST_CHECK_CLOSE(AstartDir[0], expected.momenta[0].Unit().X(), 0.01);
240  BOOST_CHECK_CLOSE(AstartDir[1], expected.momenta[0].Unit().Y(), 0.01);
241  BOOST_CHECK_CLOSE(AstartDir[2], expected.momenta[0].Unit().Z(), 0.01);
242  BOOST_CHECK_CLOSE(AendDir[0], expected.momenta[NPoints - 1].Unit().X(), 0.01);
243  BOOST_CHECK_CLOSE(AendDir[1], expected.momenta[NPoints - 1].Unit().Y(), 0.01);
244  BOOST_CHECK_CLOSE(AendDir[2], expected.momenta[NPoints - 1].Unit().Z(), 0.01);
245 
246  recob::Trajectory::Vector_t startDir, endDir;
247  std::tie(startDir, endDir) = traj.Direction();
248  BOOST_TEST_MESSAGE("Direction() start");
249  CheckVectorsClose(startDir, expected.momenta[0].Unit());
250  BOOST_TEST_MESSAGE("Direction() end");
251  CheckVectorsClose(endDir, expected.momenta[NPoints - 1].Unit());
252 
253 
254  //----------------------------------------------------------------------------
255 
256  for (size_t i = 0; i < NPoints; ++i) {
257 
258  auto const& dir = expected.momenta[i];
259  recob::Trajectory::Vector_t localDir(0., 0., dir.R());
260 
261  BOOST_TEST_MESSAGE("Test transformation to local at point #" << i);
262  auto toLocal = traj.GlobalToLocalRotationAtPoint(i);
263  CheckVectorsClose(toLocal * dir, localDir);
264 
265  BOOST_TEST_MESSAGE("Test transformation to global at point #" << i);
266  auto toGlobal = traj.LocalToGlobalRotationAtPoint(i);
267  CheckVectorsClose(toGlobal * localDir, dir);
268 
269  } // for
270 
271  //----------------------------------------------------------------------------
272  for (size_t i = 0; i < NPoints; ++i) {
273 
274  auto const& dir = expected.momenta[i];
275  recob::Trajectory::Vector_t localDir(0., 0., dir.R());
276 
277  BOOST_TEST_MESSAGE("Test legacy transformation to local at point #" << i);
278  TMatrixD TtoLocal = traj.GlobalToLocalRotationAtPoint<TMatrixD>(i);
279  auto toLocal = makeRotationMatrix(TtoLocal);
280  CheckVectorsClose(toLocal * dir, localDir);
281 
282  BOOST_TEST_MESSAGE("Test legacy transformation to global at point #" << i);
283  TMatrixD TtoGlobal = traj.LocalToGlobalRotationAtPoint<TMatrixD>(i);
284  auto toGlobal = makeRotationMatrix(TtoGlobal);
285  CheckVectorsClose(toGlobal * localDir, dir);
286 
287  } // for
288 
289  //----------------------------------------------------------------------------
290 
291 } // TestTrajectory()
292 
293 
294 //------------------------------------------------------------------------------
296 
297  BOOST_TEST_MESSAGE("Testing the default recob::Trajectory constructor");
298 
299  //
300  // Part I: initialization of trajectory inputs
301  //
302  // these are the values expected for a default-constructed trajectory
303  Expected_t expected;
304  expected.positions.clear();
305  expected.momenta.clear();
306  expected.hasMomenta = false;
307  expected.length = 0.0;
308 
309  //
310  // Part II: default constructor
311  //
312  // step II.1: create a trajectory with the default constructor
313  recob::Trajectory traj;
314 
315  for (unsigned int v = 0; v <= recob::Trajectory::MaxDumpVerbosity; ++v) {
316  std::cout << "Default-constructed trajectory dump with verbosity level "
317  << v << ":" << std::endl;
318  traj.Dump(std::cout, v, " ", " ");
319  std::cout << std::endl;
320  } // for
321 
322  // step II.2: verify that the values are as expected
323  TestTrajectory(traj, expected);
324 
325 } // TrajectoryTestDefaultConstructor()
326 
327 
328 //------------------------------------------------------------------------------
330 
331  BOOST_TEST_MESSAGE("Testing the main recob::Trajectory constructor");
332 
333  //
334  // Part I: initialization of trajectory inputs
335  //
336  const recob::Trajectory::Coord_t V2_2 = std::sqrt(0.5);
337  Expected_t expected;
338  // we describe a trajectory with uniform electric and magnetic fields aligned
339  // on z; curvature is 1 on the x/y plane.
340  expected.positions = {
341  recob::Trajectory::Point_t( -1.0, 0.0, 0.0 ),
342  recob::Trajectory::Point_t( -V2_2, +V2_2, 1.0 ),
343  recob::Trajectory::Point_t( 0.0, +1.0, 2.0 ),
344  recob::Trajectory::Point_t( +V2_2, +V2_2, 3.0 ),
345  recob::Trajectory::Point_t( +1.0, 0.0, 4.0 ),
346  recob::Trajectory::Point_t( +V2_2, -V2_2, 5.0 ),
347  recob::Trajectory::Point_t( 0.0, -1.0, 6.0 ),
348  recob::Trajectory::Point_t( -V2_2, -V2_2, 7.0 ),
349  recob::Trajectory::Point_t( -1.0, 0.0, 8.0 ),
350  recob::Trajectory::Point_t( -V2_2, +V2_2, 9.0 ),
351  recob::Trajectory::Point_t( 0.0, +1.0, 10.0 ),
352  recob::Trajectory::Point_t( +V2_2, +V2_2, 11.0 ),
353  recob::Trajectory::Point_t( +1.0, 0.0, 12.0 ),
354  recob::Trajectory::Point_t( +V2_2, -V2_2, 13.0 ),
355  recob::Trajectory::Point_t( 0.0, -1.0, 14.0 ),
356  recob::Trajectory::Point_t( -V2_2, -V2_2, 15.0 ),
357  recob::Trajectory::Point_t( -1.0, 0.0, 16.0 ),
358  recob::Trajectory::Point_t( -V2_2, +V2_2, 17.0 ),
359  recob::Trajectory::Point_t( 0.0, +1.0, 18.0 ),
360  recob::Trajectory::Point_t( +V2_2, +V2_2, 19.0 ),
361  recob::Trajectory::Point_t( +1.0, 0.0, 20.0 ),
362  recob::Trajectory::Point_t( +V2_2, -V2_2, 21.0 ),
363  recob::Trajectory::Point_t( 0.0, -1.0, 22.0 ),
364  recob::Trajectory::Point_t( -V2_2, -V2_2, 23.0 ),
365  recob::Trajectory::Point_t( -1.0, 0.0, 24.0 )
366  };
367  expected.momenta = {
368  recob::Trajectory::Vector_t( 0.0, +1.0, 1.0 ),
369  recob::Trajectory::Vector_t( +V2_2, +V2_2, 1.0 ),
370  recob::Trajectory::Vector_t( +1.0, 0.0, 1.0 ),
371  recob::Trajectory::Vector_t( +V2_2, -V2_2, 1.0 ),
372  recob::Trajectory::Vector_t( 0.0, -1.0, 1.0 ),
373  recob::Trajectory::Vector_t( -V2_2, -V2_2, 1.0 ),
374  recob::Trajectory::Vector_t( -1.0, 0.0, 1.0 ),
375  recob::Trajectory::Vector_t( -V2_2, +V2_2, 1.0 ),
376  recob::Trajectory::Vector_t( 0.0, +1.0, 1.0 ),
377  recob::Trajectory::Vector_t( +V2_2, +V2_2, 1.0 ),
378  recob::Trajectory::Vector_t( +1.0, 0.0, 1.0 ),
379  recob::Trajectory::Vector_t( +V2_2, -V2_2, 1.0 ),
380  recob::Trajectory::Vector_t( 0.0, -1.0, 1.0 ),
381  recob::Trajectory::Vector_t( -V2_2, -V2_2, 1.0 ),
382  recob::Trajectory::Vector_t( -1.0, 0.0, 1.0 ),
383  recob::Trajectory::Vector_t( -V2_2, +V2_2, 1.0 ),
384  recob::Trajectory::Vector_t( 0.0, +1.0, 1.0 ),
385  recob::Trajectory::Vector_t( +V2_2, +V2_2, 1.0 ),
386  recob::Trajectory::Vector_t( +1.0, 0.0, 1.0 ),
387  recob::Trajectory::Vector_t( +V2_2, -V2_2, 1.0 ),
388  recob::Trajectory::Vector_t( 0.0, -1.0, 1.0 ),
389  recob::Trajectory::Vector_t( -V2_2, -V2_2, 1.0 ),
390  recob::Trajectory::Vector_t( -1.0, 0.0, 1.0 ),
391  recob::Trajectory::Vector_t( -V2_2, +V2_2, 1.0 ),
392  recob::Trajectory::Vector_t( 0.0, +1.0, 1.0 )
393  };
394  expected.hasMomenta = true;
395  expected.length
396  = (expected.positions.size() - 1) * std::sqrt(3.0 - 2.0 * V2_2); //
397  expected.theta = util::pi() / 4.0;
398  expected.phi = util::pi() / 2.0;
399  expected.zenith = 0.75 * util::pi();
400  expected.azimuth = 0.0;
401 
402  //
403  // Part II: complete constructor
404  //
405  // step II.1: create a track with momentum information
406  auto positions = expected.positions;
407  auto momenta = expected.momenta;
408  recob::Trajectory traj(std::move(positions), std::move(momenta), true);
409 
410  for (unsigned int v = 0; v <= recob::Trajectory::MaxDumpVerbosity; ++v) {
411  std::cout << "Trajectory dump with verbosity level "
412  << v << ":" << std::endl;
413  traj.Dump(std::cout, v, " ", " ");
414  std::cout << std::endl;
415  } // for
416 
417  // step II.2: verify that the values are as expected
418  TestTrajectory(traj, expected);
419 
420  //
421  // Part III: complete constructor, no momentum
422  //
423 
424  // step III.1: amend the expectation for a momentumless track
425  std::transform(expected.momenta.begin(), expected.momenta.end(),
426  expected.momenta.begin(), [](auto const& v){ return v.unit(); });
427  expected.hasMomenta = false;
428 
429  // step III.2: create a track with no momentum information
430  positions = expected.positions; // copy again
431  recob::Trajectory::Momenta_t directions = expected.momenta;
432  recob::Trajectory mltraj(std::move(positions), std::move(directions), false);
433 
434  for (unsigned int v = 0; v <= recob::Trajectory::MaxDumpVerbosity; ++v) {
435  std::cout << "Momentumless trajectory dump with verbosity level "
436  << v << ":" << std::endl;
437  mltraj.Dump(std::cout, v, " ", " ");
438  std::cout << std::endl;
439  } // for
440 
441  // step III.3: verify that the values are as expected
442  TestTrajectory(mltraj, expected);
443 
444 } // TrajectoryTestMainConstructor()
445 
446 
447 //------------------------------------------------------------------------------
448 //--- registration of tests
449 //
450 // Boost needs now to know which tests we want to run.
451 // Tests are "automatically" registered, hence the BOOST_AUTO_TEST_CASE()
452 // macro name. The argument is the name of the test; each step may have a
453 // number of checks and it will fail if any of them does.
454 //
455 
456 BOOST_AUTO_TEST_CASE(TrajectoryTestDefaultConstructorTestCase) {
458 }
459 
460 BOOST_AUTO_TEST_CASE(TrajectoryTestMainConstructorTestCase) {
462 }
463 
Point_t const & Vertex() const
Returns the position of the first point of the trajectory [cm].
Definition: Trajectory.h:218
BOOST_AUTO_TEST_CASE(AllTests)
double Phi(size_t p=0) const
Azimuthal angle at a point on the trajectory, with respect to z.
Definition: Trajectory.h:337
Data product for reconstructed trajectory in space.
double EndMomentum() const
Computes and returns the modulus of momentum at the last point [GeV/c].
Definition: Trajectory.h:401
bool HasPoint(size_t i) const
Returns whether the specified trajectory point is available.
Definition: Trajectory.h:183
recob::TrackTrajectory::Positions_t positions
static constexpr Sample_t transform(Sample_t sample)
Vector_t EndDirection() const
Returns the direction of the trajectory at the last point.
Definition: Trajectory.h:303
auto const tol
Definition: SurfXYZTest.cc:16
Vector_t const & StartMomentumVector() const
Returns the momentum of the trajectory at the first point [GeV/c].
Definition: Trajectory.h:384
tracking::Coord_t Coord_t
Type used for coordinates and values in general.
Definition: Trajectory.h:70
Point_t const & End() const
Returns the position of the last point of the trajectory [cm].
Definition: Trajectory.h:226
double AzimuthAngle(size_t p=0) const
&quot;Azimuth&quot; angle of trajectory, with respect to the sky.
Definition: Trajectory.cxx:96
double ZenithAngle(size_t p=0) const
&quot;Zenith&quot; angle of trajectory, with respect to the vertical axis.
Definition: Trajectory.cxx:82
double Length(size_t startAt=0) const
Returns the approximate length of the trajectory.
Definition: Trajectory.cxx:63
Rotation_t LocalToGlobalRotationAtPoint(size_t p) const
Returns a rotation matrix bringing relative directions to global.
Definition: Trajectory.cxx:136
recob::TrackTrajectory::Rotation_t makeRotationMatrix(TMatrixD const &Tm)
double VertexMomentum() const
Computes and returns the modulus of momentum at the first point [GeV/c].
Definition: Trajectory.h:393
void CheckValue(T v, T exp, T tol, std::string tag="")
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
double StartMomentum() const
Computes and returns the modulus of momentum at the first point [GeV/c].
Definition: Trajectory.h:397
void CheckVectorsEqual(VectA const &v, VectB const &exp)
Vector_t DirectionAtPoint(size_t i) const
Computes and returns the direction of the trajectory at a point.
Definition: Trajectory.cxx:117
Rotation_t GlobalToLocalRotationAtPoint(size_t p) const
Returns a rotation matrix that brings trajectory direction along z.
Definition: Trajectory.cxx:128
recob::Trajectory::Momenta_t momenta
tracking::Vector_t Vector_t
Type for representation of momenta in 3D space.
Definition: Trajectory.h:76
T abs(T value)
size_t LastPoint() const
Returns the index of the last point in the trajectory.
Definition: Trajectory.h:175
std::pair< T, T > Extent() const
Fills the first and last point in the trajectory.
Definition: Trajectory.h:248
double MomentumAtPoint(size_t i) const
Computes and returns the modulus of the momentum at a point.
Definition: Trajectory.h:442
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
void TestTrajectory(recob::Trajectory const &traj, Expected_t const &expected)
std::pair< T, T > Direction() const
Fills the starting and ending direction of the trajectory.
Definition: Trajectory.h:472
Point_t const & LocationAtPoint(size_t i) const
Returns the position at the specified trajectory point.
Definition: Trajectory.h:236
Vector_t VertexDirection() const
Returns the direction of the trajectory at the first point.
Definition: Trajectory.h:295
recob::Trajectory::Positions_t positions
tracking::Positions_t Positions_t
Type of trajectory point list.
Definition: Trajectory.h:79
static constexpr unsigned int MaxDumpVerbosity
Largest verbosity level supported by Dump().
Definition: Trajectory.h:662
size_t NumberTrajectoryPoints() const
Returns the number of stored trajectory points.
Definition: Trajectory.h:156
recob::TrackTrajectory::Momenta_t momenta
tracking::Momenta_t Momenta_t
Type of momentum list.
Definition: Trajectory.h:82
tuple dir
Definition: dropbox.py:28
Vector_t const & MomentumVectorAtPoint(size_t i) const
Returns the momentum vector at a point.
Definition: Trajectory.h:457
Vector_t StartDirection() const
Returns the direction of the trajectory at the first point.
Definition: Trajectory.h:299
A trajectory in space reconstructed from hits.
Definition: Trajectory.h:67
bool HasMomentum() const
Returns whether information about the momentum is available.
Definition: Trajectory.h:425
tracking::Point_t Point_t
Type for representation of position in physical 3D space.
Definition: Trajectory.h:73
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
void TrajectoryTestDefaultConstructor()
Vector_t const & VertexMomentumVector() const
Returns the momentum of the trajectory at the first point [GeV/c].
Definition: Trajectory.h:380
double Theta(size_t p=0) const
Trajectory angle at point, with respect to positive z direction.
Definition: Trajectory.h:320
void Dump(Stream &&out, unsigned int verbosity, std::string indent, std::string indentFirst) const
Prints trajectory content into a stream.
Collection of Physical constants used in LArSoft.
Point_t const & Start() const
Returns the position of the first point of the trajectory [cm].
Definition: Trajectory.h:222
tracking::Rotation_t Rotation_t
Type for representation of space rotations.
Definition: Trajectory.h:97
BEGIN_PROLOG could also be cout
void TrajectoryTestMainConstructor()
Vector_t const & EndMomentumVector() const
Returns the momentum of the trajectory at the last point [GeV/c].
Definition: Trajectory.h:388
void CheckVectorsClose(VectA const &v, VectB const &exp, double tol=0.01)