LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
trkf::PropYZLine Class Reference

#include "PropYZLine.h"

Inheritance diagram for trkf::PropYZLine:
trkf::Propagator

Public Types

enum  PropDirection { FORWARD, BACKWARD, UNKNOWN }
 Propagation direction enum. More...
 

Public Member Functions

 PropYZLine (double tcut, bool doDedx)
 Constructor. More...
 
virtual ~PropYZLine ()
 Destructor. More...
 
Propagatorclone () const
 Clone method. More...
 
boost::optional< double > short_vec_prop (KTrack &trk, const std::shared_ptr< const Surface > &surf, Propagator::PropDirection dir, bool doDedx, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const
 Propagate without error. More...
 
virtual boost::optional< double > origin_vec_prop (KTrack &trk, const std::shared_ptr< const Surface > &porient, TrackMatrix *prop_matrix=0) const
 Propagate without error to surface whose origin parameters coincide with track position. More...
 
double getTcut () const
 
bool getDoDedx () const
 
const std::shared_ptr< const Interactor > & getInteractor () const
 
boost::optional< double > vec_prop (KTrack &trk, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const
 Propagate without error (long distance). More...
 
boost::optional< double > lin_prop (KTrack &trk, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, KTrack *ref=0, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const
 Linearized propagate without error. More...
 
boost::optional< double > err_prop (KETrack &tre, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, KTrack *ref=0, TrackMatrix *prop_matrix=0) const
 Propagate with error, but without noise. More...
 
boost::optional< double > noise_prop (KETrack &tre, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, KTrack *ref=0) const
 Propagate with error and noise. More...
 
boost::optional< double > dedx_prop (double pinv, double mass, double s, double *deriv=0) const
 Method to calculate updated momentum due to dE/dx. More...
 

Private Member Functions

bool transformYZLine (double phi1, double phi2, TrackVector &vec, Surface::TrackDirection &dir, TrackMatrix *prop_matrix) const
 Transform yz line -> yz line. More...
 
bool transformYZPlane (double phi1, double phi2, TrackVector &vec, Surface::TrackDirection &dir, TrackMatrix *prop_matrix) const
 Transform yz plane -> yz line. More...
 
bool transformXYZPlane (double theta1, double phi1, double phi2, TrackVector &vec, Surface::TrackDirection &dir, TrackMatrix *prop_matrix) const
 Transform xyz plane -> yz line. More...
 

Detailed Description

Definition at line 20 of file PropYZLine.h.

Member Enumeration Documentation

Propagation direction enum.

Enumerator
FORWARD 
BACKWARD 
UNKNOWN 

Definition at line 92 of file Propagator.h.

Constructor & Destructor Documentation

trkf::PropYZLine::PropYZLine ( double  tcut,
bool  doDedx 
)

Constructor.

Constructor.

Arguments.

tcut - Delta ray energy cutoff for calculating dE/dx. doDedx - dE/dx enable flag.

Definition at line 28 of file PropYZLine.cxx.

Referenced by clone().

28  :
29  Propagator(tcut, doDedx, (tcut >= 0. ?
30  std::shared_ptr<const Interactor>(new InteractGeneral(tcut)) :
31  std::shared_ptr<const Interactor>()))
32  {}
Propagator(double tcut, bool doDedx, const std::shared_ptr< const Interactor > &interactor)
Constructor.
Definition: Propagator.cxx:25
trkf::PropYZLine::~PropYZLine ( )
virtual

Destructor.

Definition at line 35 of file PropYZLine.cxx.

36  {}

Member Function Documentation

Propagator* trkf::PropYZLine::clone ( ) const
inlinevirtual

Clone method.

Implements trkf::Propagator.

Definition at line 33 of file PropYZLine.h.

References dir, origin_vec_prop(), PropYZLine(), short_vec_prop(), transformXYZPlane(), transformYZLine(), and transformYZPlane().

33 {return new PropYZLine(*this);}
PropYZLine(double tcut, bool doDedx)
Constructor.
Definition: PropYZLine.cxx:28
boost::optional< double > trkf::Propagator::dedx_prop ( double  pinv,
double  mass,
double  s,
double *  deriv = 0 
) const
inherited

Method to calculate updated momentum due to dE/dx.

Method to calculate updated momentum due to dE/dx.

Arguments:

pinv - Initial inverse momentum (units c/GeV). mass - Particle mass (GeV/c^2). s - Path distance. deriv - Pointer to store derivative d(pinv2)/d(pinv1) if nonzero.

Returns: Final inverse momentum (pinv2) + success flag.

Failure is returned in case of range out.

Inverse momentum can be signed (q/p). Returned inverse momentum has the same sign as the input.

In this method, we are solving the differential equation in terms of energy.

dE/dx = -f(E)

where f(E) is the stopping power returned by method LArProperties::Eloss.

We expect that this method will be called exclusively for short distance propagation. The differential equation is solved using the midpoint method using a single step, which requires two evaluations of f(E).

dE = -s*f(E1) E2 = E1 - s*f(E1 + 0.5*dE)

The derivative is calculated assuming E2 = E1 + constant, giving

d(pinv2)/d(pinv1) = pinv2^3 E2 / (pinv1^3 E1).

Definition at line 459 of file Propagator.cxx.

References trkf::Propagator::fTcut.

Referenced by trkf::Propagator::getInteractor(), trkf::PropXYZPlane::short_vec_prop(), short_vec_prop(), and trkf::PropYZPlane::short_vec_prop().

461  {
462  // For infinite initial momentum, return with success status,
463  // still infinite momentum.
464 
465  if(pinv == 0.)
466  return boost::optional<double>(true, 0.);
467 
468  // Set the default return value to be uninitialized with value 0.
469 
470  boost::optional<double> result(false, 0.);
471 
472  // Get LAr service.
473 
474  auto const * detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
475 
476  // Calculate final energy.
477 
478  double p1 = 1./std::abs(pinv);
479  double e1 = std::hypot(p1, mass);
480  double de = -0.001 * s * detprop->Eloss(p1, mass, fTcut);
481  double emid = e1 + 0.5 * de;
482  if(emid > mass) {
483  double pmid = std::sqrt(emid*emid - mass*mass);
484  double e2 = e1 - 0.001 * s * detprop->Eloss(pmid, mass, fTcut);
485  if(e2 > mass) {
486  double p2 = std::sqrt(e2*e2 - mass*mass);
487  double pinv2 = 1./p2;
488  if(pinv < 0.)
489  pinv2 = -pinv2;
490 
491  // Calculation was successful, update result.
492 
493  result = boost::optional<double>(true, pinv2);
494 
495  // Also calculate derivative, if requested.
496 
497  if(deriv != 0)
498  *deriv = pinv2*pinv2*pinv2 * e2 / (pinv*pinv*pinv * e1);
499  }
500  }
501 
502  // Done.
503 
504  return result;
505  }
Float_t s
Definition: plot.C:23
double fTcut
Maximum delta ray energy for dE/dx.
Definition: Propagator.h:163
boost::optional< double > trkf::Propagator::err_prop ( KETrack tre,
const std::shared_ptr< const Surface > &  psurf,
PropDirection  dir,
bool  doDedx,
KTrack ref = 0,
TrackMatrix prop_matrix = 0 
) const
inherited

Propagate with error, but without noise.

Propagate with error, but without noise (i.e. reversibly).

Arguments:

tre - Track to propagate. psurf - Destination surface. dir - Propagation direction (FORWARD, BACKWARD, or UNKNOWN). doDedx - dE/dx enable/disable flag. ref - Reference track (for linearized propagation). Can be null. prop_matrix - Return propagation matrix if not null.

Returned value: propagation distance + success flag.

Definition at line 354 of file Propagator.cxx.

References trkf::KETrack::getError(), trkf::Propagator::lin_prop(), and trkf::KETrack::setError().

Referenced by trkf::KGTrack::fillTrack(), trkf::KalmanFilterAlg::fitMomentumMS(), trkf::Propagator::getInteractor(), and trkf::KHit< N >::predict().

360  {
361  // Propagate without error, get propagation matrix.
362 
363  TrackMatrix prop_temp;
364  if(prop_matrix == 0)
365  prop_matrix = &prop_temp;
366  boost::optional<double> result = lin_prop(tre, psurf, dir, doDedx, ref, prop_matrix, 0);
367 
368  // If propagation succeeded, update track error matrix.
369 
370  if(!!result) {
371  TrackMatrix temp = prod(tre.getError(), trans(*prop_matrix));
372  TrackMatrix temp2 = prod(*prop_matrix, temp);
373  TrackError newerr = ublas::symmetric_adaptor<TrackMatrix>(temp2);
374  tre.setError(newerr);
375  }
376 
377  // Done.
378 
379  return result;
380  }
KSymMatrix< 5 >::type TrackError
Track error matrix, dimension 5x5.
TDirectory * dir
Definition: macro.C:5
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
boost::optional< double > lin_prop(KTrack &trk, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, KTrack *ref=0, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const
Linearized propagate without error.
Definition: Propagator.cxx:258
bool trkf::Propagator::getDoDedx ( ) const
inlineinherited
const std::shared_ptr<const Interactor>& trkf::Propagator::getInteractor ( ) const
inlineinherited
double trkf::Propagator::getTcut ( ) const
inlineinherited

Definition at line 102 of file Propagator.h.

References trkf::Propagator::fTcut.

102 {return fTcut;}
double fTcut
Maximum delta ray energy for dE/dx.
Definition: Propagator.h:163
boost::optional< double > trkf::Propagator::lin_prop ( KTrack trk,
const std::shared_ptr< const Surface > &  psurf,
PropDirection  dir,
bool  doDedx,
KTrack ref = 0,
TrackMatrix prop_matrix = 0,
TrackError noise_matrix = 0 
) const
inherited

Linearized propagate without error.

Linearized propagate without error.

Arguments:

trk - Track to propagate. psurf - Destination surface. dir - Propagation direction (FORWARD, BACKWARD, or UNKNOWN). doDedx - dE/dx enable/disable flag. ref - Reference track (for linearized propagation). Can be null. prop_matrix - Return propagation matrix if not null. noise_matrix - Return noise matrix if not null.

Returned value: Propagation distance & success flag.

If the reference track is null, this method simply calls vec_prop.

Definition at line 258 of file Propagator.cxx.

References trkf::KTrack::getDirection(), trkf::KTrack::getSurface(), trkf::KTrack::getVector(), trkf::KTrack::isValid(), trkf::KTrack::setDirection(), trkf::KTrack::setSurface(), trkf::KTrack::setVector(), and trkf::Propagator::vec_prop().

Referenced by trkf::Propagator::err_prop(), trkf::Propagator::getInteractor(), and trkf::Propagator::noise_prop().

265  {
266  // Default result.
267 
268  boost::optional<double> result;
269 
270  if(ref == 0)
271  result = vec_prop(trk, psurf, dir, doDedx, prop_matrix, noise_matrix);
272  else {
273 
274  // A reference track has been provided.
275 
276  // It is an error (throw exception) if the reference track and
277  // the track to be propagted are not on the same surface.
278 
279  if(!trk.getSurface()->isEqual(*(ref->getSurface())))
280  throw cet::exception("Propagator") <<
281  "Input track and reference track not on same surface.\n";
282 
283  // Remember the starting track and reference track.
284 
285  KTrack trk0(trk);
286  KTrack ref0(*ref);
287 
288  // Propagate the reference track. Make sure we calculate the
289  // propagation matrix.
290 
291  TrackMatrix prop_temp;
292  if(prop_matrix == 0)
293  prop_matrix = &prop_temp;
294 
295  // Do the propgation. The returned result will be the result of
296  // this propagatrion.
297 
298  result = vec_prop(*ref, psurf, dir, doDedx, prop_matrix, noise_matrix);
299  if(!!result) {
300 
301  // Propagation of reference track succeeded. Update the track
302  // state vector and surface of the track to be propagated.
303 
304  TrackVector diff = trk.getSurface()->getDiff(trk.getVector(), ref0.getVector());
305  TrackVector newvec = ref->getVector() + prod(*prop_matrix, diff);
306 
307  // Store updated state vector and surface.
308 
309  trk.setVector(newvec);
310  trk.setSurface(psurf);
311  trk.setDirection(ref->getDirection());
312  if (!trk.getSurface()->isEqual(*(ref->getSurface())))
313  throw cet::exception("Propagator") << __func__ << ": surface mismatch";
314 
315  // Final validity check. In case of failure, restore the track
316  // and reference track to their starting values.
317 
318  if(!trk.isValid()) {
319  result = boost::optional<double>(false, 0.);
320  trk = trk0;
321  *ref = ref0;
322  }
323  }
324  else {
325 
326  // Propagation failed.
327  // Restore the reference track to its starting value, so that we ensure
328  // the reference track and the actual track remain on the same surface.
329 
330  trk = trk0;
331  *ref = ref0;
332  }
333  }
334 
335  // Done.
336 
337  return result;
338  }
KVector< 5 >::type TrackVector
Track state vector, dimension 5.
boost::optional< double > vec_prop(KTrack &trk, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const
Propagate without error (long distance).
Definition: Propagator.cxx:52
TDirectory * dir
Definition: macro.C:5
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
boost::optional< double > trkf::Propagator::noise_prop ( KETrack tre,
const std::shared_ptr< const Surface > &  psurf,
PropDirection  dir,
bool  doDedx,
KTrack ref = 0 
) const
inherited

Propagate with error and noise.

Propagate with error and noise.

Arguments:

tre - Track to propagate. psurf - Destination surface. dir - Propagation direction (FORWARD, BACKWARD, or UNKNOWN). doDedx - dE/dx enable/disable flag. ref - Reference track (for linearized propagation). Can be null.

Returned value: propagation distance + success flag.

Definition at line 394 of file Propagator.cxx.

References trkf::KETrack::getError(), trkf::Propagator::lin_prop(), and trkf::KETrack::setError().

Referenced by trkf::KalmanFilterAlg::buildTrack(), trkf::KalmanFilterAlg::extendTrack(), trkf::KalmanFilterAlg::fitMomentumMS(), trkf::Propagator::getInteractor(), trkf::KalmanFilterAlg::smoothTrack(), and trkf::KalmanFilterAlg::updateMomentum().

399  {
400  // Propagate without error, get propagation matrix and noise matrix.
401 
402  TrackMatrix prop_matrix;
403  TrackError noise_matrix;
404  boost::optional<double> result = lin_prop(tre, psurf, dir, doDedx, ref,
405  &prop_matrix, &noise_matrix);
406 
407  // If propagation succeeded, update track error matrix.
408 
409  if(!!result) {
410  TrackMatrix temp = prod(tre.getError(), trans(prop_matrix));
411  TrackMatrix temp2 = prod(prop_matrix, temp);
412  TrackError newerr = ublas::symmetric_adaptor<TrackMatrix>(temp2);
413  newerr += noise_matrix;
414  tre.setError(newerr);
415  }
416 
417  // Done.
418 
419  return result;
420  }
KSymMatrix< 5 >::type TrackError
Track error matrix, dimension 5x5.
TDirectory * dir
Definition: macro.C:5
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
boost::optional< double > lin_prop(KTrack &trk, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, KTrack *ref=0, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const
Linearized propagate without error.
Definition: Propagator.cxx:258
boost::optional< double > trkf::PropYZLine::origin_vec_prop ( KTrack trk,
const std::shared_ptr< const Surface > &  porient,
TrackMatrix prop_matrix = 0 
) const
virtual

Propagate without error to surface whose origin parameters coincide with track position.

Propagate without error to dynamically generated origin surface. Optionally return propagation matrix.

Arguments:

trk - Track to propagate. porient - Orientation surface. prop_matrix - Pointer to optional propagation matrix.

Returned value: propagation distance + success flag.

Propagation distance is always zero after successful propagation.

Implements trkf::Propagator.

Definition at line 272 of file PropYZLine.cxx.

References dir, trkf::KTrack::getDirection(), trkf::KTrack::getPosition(), trkf::KTrack::getSurface(), trkf::KTrack::getVector(), trkf::KTrack::isValid(), trkf::SurfYZLine::phi(), trkf::KTrack::setDirection(), trkf::KTrack::setSurface(), trkf::KTrack::setVector(), transformXYZPlane(), transformYZLine(), and transformYZPlane().

Referenced by clone(), trkf::PropAny::origin_vec_prop(), and short_vec_prop().

275  {
276  // Set the default return value to be unitialized with value 0.
277 
278  boost::optional<double> result(false, 0.);
279 
280  // Remember starting track.
281 
282  KTrack trk0(trk);
283 
284  // Get initial track parameters and direction.
285  // Note the initial track can be on any type of surface.
286 
287  TrackVector vec = trk.getVector(); // Modifiable copy.
288  if(vec.size() != 5)
289  throw cet::exception("PropYZPlane")
290  << "Track state vector has wrong size" << vec.size() << "\n";
291  Surface::TrackDirection dir = trk.getDirection();
292 
293  // Get track position.
294 
295  double xyz[3];
296  trk.getPosition(xyz);
297  double x02 = xyz[0];
298  double y02 = xyz[1];
299  double z02 = xyz[2];
300 
301  // Generate the origin surface, which will be the destination surface.
302  // Return failure if orientation surface is the wrong type.
303 
304  const SurfYZLine* orient = dynamic_cast<const SurfYZLine*>(&*porient);
305  if(orient == 0)
306  return result;
307  double phi2 = orient->phi();
308  std::shared_ptr<const Surface> porigin(new SurfYZLine(x02, y02, z02, phi2));
309 
310  // Test initial surface types.
311 
312  if(const SurfYZLine* from = dynamic_cast<const SurfYZLine*>(&*trk.getSurface())) {
313 
314  // Initial surface is SurfYZLine.
315  // Get surface paramters.
316 
317  double phi1 = from->phi();
318 
319  // Transform track to origin surface.
320 
321  bool ok = transformYZLine(phi1, phi2, vec, dir, prop_matrix);
322  result = boost::optional<double>(ok, 0.);
323  if(!ok)
324  return result;
325  }
326  else if(const SurfYZPlane* from = dynamic_cast<const SurfYZPlane*>(&*trk.getSurface())) {
327 
328  // Initial surface is SurfYZPlane.
329  // Get surface paramters.
330 
331  double phi1 = from->phi();
332 
333  // Transform track to origin surface.
334 
335  bool ok = transformYZPlane(phi1, phi2, vec, dir, prop_matrix);
336  result = boost::optional<double>(ok, 0.);
337  if(!ok)
338  return result;
339  }
340  else if(const SurfXYZPlane* from = dynamic_cast<const SurfXYZPlane*>(&*trk.getSurface())) {
341 
342  // Initial surface is SurfXYZPlane.
343  // Get surface paramters.
344 
345  double theta1 = from->theta();
346  double phi1 = from->phi();
347 
348  // Transform track to origin surface.
349 
350  bool ok = transformXYZPlane(theta1, phi1, phi2, vec, dir, prop_matrix);
351  result = boost::optional<double>(ok, 0.);
352  if(!ok)
353  return result;
354  }
355 
356  // Update track.
357 
358  trk.setSurface(porigin);
359  trk.setVector(vec);
360  trk.setDirection(dir);
361 
362  // Final validity check.
363 
364  if(!trk.isValid()) {
365  trk = trk0;
366  result = boost::optional<double>(false, 0.);
367  }
368 
369  // Done.
370 
371  return result;
372  }
TrackDirection
Track direction enum.
Definition: Surface.h:56
bool transformXYZPlane(double theta1, double phi1, double phi2, TrackVector &vec, Surface::TrackDirection &dir, TrackMatrix *prop_matrix) const
Transform xyz plane -> yz line.
Definition: PropYZLine.cxx:737
bool transformYZLine(double phi1, double phi2, TrackVector &vec, Surface::TrackDirection &dir, TrackMatrix *prop_matrix) const
Transform yz line -> yz line.
Definition: PropYZLine.cxx:376
bool transformYZPlane(double phi1, double phi2, TrackVector &vec, Surface::TrackDirection &dir, TrackMatrix *prop_matrix) const
Transform yz plane -> yz line.
Definition: PropYZLine.cxx:564
KVector< 5 >::type TrackVector
Track state vector, dimension 5.
TDirectory * dir
Definition: macro.C:5
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
boost::optional< double > trkf::PropYZLine::short_vec_prop ( KTrack trk,
const std::shared_ptr< const Surface > &  psurf,
Propagator::PropDirection  dir,
bool  doDedx,
TrackMatrix prop_matrix = 0,
TrackError noise_matrix = 0 
) const
virtual

Propagate without error.

Propagate without error. Optionally return propagation matrix and noise matrix.

Arguments:

trk - Track to propagate. psurf - Destination surface. dir - Propagation direction (FORWARD, BACKWARD, or UNKNOWN). doDedx - dE/dx enable/disable flag. prop_matrix - Pointer to optional propagation matrix. noise_matrix - Pointer to optional noise matrix.

Returned value: propagation distance + success flag.

Implements trkf::Propagator.

Definition at line 53 of file PropYZLine.cxx.

References trkf::Propagator::BACKWARD, trkf::Propagator::dedx_prop(), trkf::Propagator::FORWARD, trkf::Propagator::getDoDedx(), trkf::Propagator::getInteractor(), trkf::KTrack::getPosition(), trkf::KTrack::getVector(), trkf::KTrack::Mass(), origin_vec_prop(), trkf::SurfYZLine::phi(), s, trkf::KTrack::setSurface(), trkf::KTrack::setVector(), trkf::Propagator::UNKNOWN, trkf::SurfYZLine::x0(), trkf::SurfYZLine::y0(), and trkf::SurfYZLine::z0().

Referenced by clone(), and trkf::PropAny::short_vec_prop().

59  {
60  // Set the default return value to be unitialized with value 0.
61 
62  boost::optional<double> result(false, 0.);
63 
64  // Get destination surface and surface parameters.
65  // Return failure if wrong surface type.
66 
67  const SurfYZLine* to = dynamic_cast<const SurfYZLine*>(&*psurf);
68  if(to == 0)
69  return result;
70  double x02 = to->x0();
71  double y02 = to->y0();
72  double z02 = to->z0();
73  double phi2 = to->phi();
74 
75  // Remember starting track.
76 
77  KTrack trk0(trk);
78 
79  // Get track position.
80 
81  double xyz[3];
82  trk.getPosition(xyz);
83  double x01 = xyz[0];
84  double y01 = xyz[1];
85  double z01 = xyz[2];
86 
87  // Propagate to origin surface.
88 
89  TrackMatrix local_prop_matrix;
90  TrackMatrix* plocal_prop_matrix = (prop_matrix==0 ? 0 : &local_prop_matrix);
91  boost::optional<double> result1 = origin_vec_prop(trk, psurf, plocal_prop_matrix);
92  if(!result1)
93  return result1;
94 
95  // Get the intermediate track state vector and track parameters.
96 
97  const TrackVector& vec = trk.getVector();
98  if(vec.size() != 5)
99  throw cet::exception("PropYZLine")
100  << "Track state vector has wrong size" << vec.size() << "\n";
101  double r1 = vec(0);
102  double v1 = vec(1);
103  double phid1 = vec(2);
104  double eta1 = vec(3);
105  double pinv = vec(4);
106 
107  // Calculate transcendental functions.
108 
109  double sinphid1 = std::sin(phid1);
110  double cosphid1 = std::cos(phid1);
111  double sh1 = std::sinh(eta1);
112  double ch1 = std::cosh(eta1);
113  double sinphi2 = std::sin(phi2);
114  double cosphi2 = std::cos(phi2);
115 
116  // Calculate the initial position in the intermediate coordinate system.
117 
118  double u1 = -r1 * sinphid1;
119  double w1 = r1 * cosphid1;
120 
121  // Calculate initial position in the destination coordinate system.
122 
123  double u2 = x01 - x02 + u1;
124  double v2 = (y01 - y02) * cosphi2 + (z01 - z02) * sinphi2 + v1;
125  double w2 = -(y01 - y02) * sinphi2 + (z01 - z02) * cosphi2 + w1;
126 
127  // Calculate the impact parameter in the destination coordinate system.
128 
129  double r2 = w2 * cosphid1 - u2 * sinphid1;
130 
131  // Calculate the perpendicular propagation distance.
132 
133  double d2 = -(w2 * sinphid1 + u2 * cosphid1);
134 
135  // Calculate the final position in the destination coordinate system.
136 
137  //double u2p = -r2 * sinphid1;
138  double v2p = v2 + d2 * sh1;
139  //double w2p = r2 * cosphid1;
140 
141  // Calculate the signed propagation distance.
142 
143  double s = d2 * ch1;
144 
145  // Check if propagation was in the right direction.
146  // (Compare sign of s with requested direction).
147 
148  bool sok = (dir == Propagator::UNKNOWN ||
149  (dir == Propagator::FORWARD && s >= 0.) ||
150  (dir == Propagator::BACKWARD && s <= 0.));
151 
152  // If wrong direction, return failure without updating the track
153  // or propagation matrix.
154 
155  if(!sok) {
156  trk = trk0;
157  return result;
158  }
159 
160  // Find final momentum.
161 
162  double deriv = 1.;
163  boost::optional<double> pinv2(true, pinv);
164  if(getDoDedx() && doDedx && s != 0.) {
165  double* pderiv = (prop_matrix != 0 ? &deriv : 0);
166  pinv2 = dedx_prop(pinv, trk.Mass(), s, pderiv);
167  }
168 
169  // Return failure in case of range out.
170 
171  if(!pinv2) {
172  trk = trk0;
173  return result;
174  }
175 
176  // Update default result to success and store propagation distance.
177 
178  result = boost::optional<double>(true, s);
179 
180  // Update propagation matrix (if requested).
181 
182  if(prop_matrix != 0) {
183  TrackMatrix pm;
184  pm.resize(vec.size(), vec.size(), false);
185 
186  // Calculate partial derivatives.
187 
188  pm(0,0) = 1.; // dr2/dr1
189  pm(1,0) = 0.; // dv2/dr1
190  pm(2,0) = 0.; // d(phi2)/dr1
191  pm(3,0) = 0.; // d(eta2)/dr1
192  pm(4,0) = 0.; // d(pinv2)/dr1
193 
194  pm(0,1) = 0.; // dr2/dv1
195  pm(1,1) = 1.; // dv2/dv1
196  pm(2,1) = 0.; // d(phi2)/dv1
197  pm(3,1) = 0.; // d(eta2)/dv1
198  pm(4,1) = 0.; // d(pinv2)/dv1
199 
200  pm(0,2) = d2; // dr2/d(phi1);
201  pm(1,2) = -r2*sh1; // dv2/d(phi1);
202  pm(2,2) = 1.; // d(phi2)/d(phi1);
203  pm(3,2) = 0.; // d(eta2)/d(phi1);
204  pm(4,2) = 0.; // d(pinv2)/d(phi1);
205 
206  pm(0,3) = 0.; // dr2/d(eta1);
207  pm(1,3) = d2*ch1; // dv2/d(eta1);
208  pm(2,3) = 0.; // d(phi2)/d(eta1);
209  pm(3,3) = 1.; // d(eta2)/d(eta1);
210  pm(4,3) = 0.; // d(pinv2)/d(eta1);
211 
212  pm(0,4) = 0.; // dr2/d(pinv1);
213  pm(1,4) = 0.; // dv2/d(pinv1);
214  pm(2,4) = 0.; // d(phi2)/d(pinv1);
215  pm(3,4) = 0.; // d(eta2)/d(pinv1);
216  pm(4,4) = deriv; // d(pinv2)/d(pinv1);
217 
218  // Compose the final propagation matrix from zero-distance propagation and
219  // parallel surface propagation.
220 
221  *prop_matrix = prod(pm, *plocal_prop_matrix);
222  }
223 
224  // Update noise matrix (if requested).
225 
226  if(noise_matrix != 0) {
227  noise_matrix->resize(vec.size(), vec.size(), false);
228  if(getInteractor().get() != 0) {
229  bool ok = getInteractor()->noise(trk, s, *noise_matrix);
230  if(!ok) {
231  trk = trk0;
232  return boost::optional<double>(false, 0.);
233  }
234  }
235  else
236  noise_matrix->clear();
237  }
238 
239  // Construct track vector at destination surface.
240 
241  TrackVector vec2(vec.size());
242  vec2(0) = r2;
243  vec2(1) = v2p;
244  vec2(2) = phid1;
245  vec2(3) = eta1;
246  vec2(4) = *pinv2;
247 
248  // Update track.
249 
250  trk.setSurface(psurf);
251  trk.setVector(vec2);
252 
253  // Done.
254 
255  return result;
256  }
Float_t s
Definition: plot.C:23
virtual boost::optional< double > origin_vec_prop(KTrack &trk, const std::shared_ptr< const Surface > &porient, TrackMatrix *prop_matrix=0) const
Propagate without error to surface whose origin parameters coincide with track position.
Definition: PropYZLine.cxx:272
boost::optional< double > dedx_prop(double pinv, double mass, double s, double *deriv=0) const
Method to calculate updated momentum due to dE/dx.
Definition: Propagator.cxx:459
const std::shared_ptr< const Interactor > & getInteractor() const
Definition: Propagator.h:104
KVector< 5 >::type TrackVector
Track state vector, dimension 5.
bool getDoDedx() const
Definition: Propagator.h:103
TDirectory * dir
Definition: macro.C:5
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
bool trkf::PropYZLine::transformXYZPlane ( double  theta1,
double  phi1,
double  phi2,
TrackVector vec,
Surface::TrackDirection dir,
TrackMatrix prop_matrix 
) const
private

Transform xyz plane -> yz line.

Definition at line 737 of file PropYZLine.cxx.

References trkf::Surface::BACKWARD, and trkf::Surface::FORWARD.

Referenced by clone(), and origin_vec_prop().

741  {
742  // Calculate surface transcendental functions.
743 
744  double sinth1 = std::sin(theta1);
745  double costh1 = std::cos(theta1);
746 
747  double sindphi = std::sin(phi2 - phi1);
748  double cosdphi = std::cos(phi2 - phi1);
749 
750  // Get the initial track parameters.
751 
752  double dudw1 = vec(2);
753  double dvdw1 = vec(3);
754 
755  // Make sure initial track has a valid direction.
756 
757  double dirf = 1.;
758  if(dir == Surface::BACKWARD)
759  dirf = -1.;
760  else if(dir != Surface::FORWARD)
761  return false;
762 
763  // Calculate elements of rotation matrix from initial coordinate
764  // system to destination coordinte system.
765 
766  double ruu = costh1;
767  double ruw = sinth1;
768 
769  double rvu = -sinth1*sindphi;
770  double rvv = cosdphi;
771  double rvw = costh1*sindphi;
772 
773  double rwu = -sinth1*cosdphi;
774  double rwv = -sindphi;
775  double rww = costh1*cosdphi;
776 
777  // Calculate direction in the starting coordinate system.
778 
779  double dw1 = dirf / std::sqrt(1. + dudw1*dudw1 + dvdw1*dvdw1);
780  double du1 = dudw1 * dw1;
781  double dv1 = dvdw1 * dw1;
782 
783  // Rotate direction vector into destination coordinate system.
784 
785  double du2 = ruu*du1 + ruw*dw1;
786  double dv2 = rvu*du1 + rvv*dv1 + rvw*dw1;
787  double dw2 = rwu*du1 + rwv*dv1 + rww*dw1;
788  double duw2 = std::hypot(du2, dw2);
789 
790  // Calculate final direction track parameters.
791 
792  double phid2 = atan2(dw2, du2);
793  double eta2 = std::asinh(dv2 / duw2);
794 
795  // Update propagation matrix (if requested).
796 
797  if(prop_matrix != 0) {
798  TrackMatrix& pm = *prop_matrix;
799  pm.resize(vec.size(), vec.size(), false);
800 
801  // Calculate partial derivatives.
802 
803  // Partials of initial positions and directions wrt initial t.p.'s.
804 
805  double ddu1ddudw1 = (1. + dvdw1*dvdw1) * dw1*dw1*dw1;
806  double ddu1ddvdw1 = -dudw1 * dvdw1 * dw1*dw1*dw1;
807 
808  double ddv1ddudw1 = -dudw1 * dvdw1 * dw1*dw1*dw1;
809  double ddv1ddvdw1 = (1. + dudw1*dudw1) * dw1*dw1*dw1;
810 
811  double ddw1ddudw1 = -dudw1 * dw1*dw1*dw1;
812  double ddw1ddvdw1 = -dvdw1 * dw1*dw1*dw1;
813 
814  // Rotate partials to destination coordinate system.
815 
816  double du2du1 = ruu;
817  double dv2du1 = rvu;
818  double dw2du1 = rwu;
819 
820  double dv2dv1 = rvv;
821  double dw2dv1 = rwv;
822 
823  double ddu2ddudw1 = ruu*ddu1ddudw1 + ruw*ddw1ddudw1;
824  double ddv2ddudw1 = rvu*ddu1ddudw1 + rvv*ddv1ddudw1 + rvw*ddw1ddudw1;
825  double ddw2ddudw1 = rwu*ddu1ddudw1 + rwv*ddv1ddudw1 + rww*ddw1ddudw1;
826 
827  double ddu2ddvdw1 = ruu*ddu1ddvdw1 + ruw*ddw1ddvdw1;
828  double ddv2ddvdw1 = rvu*ddu1ddvdw1 + rvv*ddv1ddvdw1 + rvw*ddw1ddvdw1;
829  double ddw2ddvdw1 = rwu*ddu1ddvdw1 + rwv*ddv1ddvdw1 + rww*ddw1ddvdw1;
830 
831  // Partials of final t.p. wrt final position and direction.
832 
833  double dr2du2 = -dw2/duw2;
834  double dr2dw2 = du2/duw2;
835 
836  double dphi2ddu2 = -dw2/(duw2*duw2);
837  double dphi2ddw2 = du2/(duw2*duw2);
838 
839  double deta2ddv2 = 1./(duw2*duw2);
840 
841  // Partials of final t.p. wrt initial t.p.
842 
843  double dr2du1 = dr2du2*du2du1 + dr2dw2*dw2du1;
844  double dr2dv1 = dr2dw2*dw2dv1;
845 
846  double dphi2ddudw1 = dphi2ddu2*ddu2ddudw1 + dphi2ddw2*ddw2ddudw1;
847  double dphi2ddvdw1 = dphi2ddu2*ddu2ddvdw1 + dphi2ddw2*ddw2ddvdw1;
848 
849  double deta2ddudw1 = deta2ddv2*ddv2ddudw1;
850  double deta2ddvdw1 = deta2ddv2*ddv2ddvdw1;
851 
852  // We still need to calculate the corretion due to the dependence of the
853  // propagation distance on the initial track parameters. This correction is
854  // needed even though the actual propagation distance is zero.
855 
856  // This correction only effects the v track parameter, since the v parameter
857  // the only parameter that actually dependes on the propagation distance.
858 
859  // Partials of propagation distance wrt position and direction in the destination
860  // coordinate system.
861 
862  double dsdu2 = -du2/(duw2*duw2);
863  double dsdw2 = -dw2/(duw2*duw2);
864 
865  // Partials of propagation distance wrt initial t.p.
866 
867  double dsdu1 = dsdu2*du2du1 + dsdw2*dw2du1;
868  double dsdv1 = dsdw2*dw2dv1;
869 
870  // Calculate correction to v parameter partials wrt initial t.p. due to path length.
871 
872  dv2du1 += dv2*dsdu1;
873  dv2dv1 += dv2*dsdv1;
874 
875  // Fill matrix.
876 
877  pm(0,0) = dr2du1; // dr2/du1
878  pm(1,0) = dv2du1; // dv2/du1
879  pm(2,0) = 0.; // d(phi2)/du1
880  pm(3,0) = 0.; // d(eta2)/du1
881  pm(4,0) = 0.; // d(pinv2)/du1
882 
883  pm(0,1) = dr2dv1; // dr2/dv1
884  pm(1,1) = dv2dv1; // dv2/dv1
885  pm(2,1) = 0.; // d(phi2)/dv1
886  pm(3,1) = 0.; // d(eta2)/dv1
887  pm(4,1) = 0.; // d(pinv2)/dv1
888 
889  pm(0,2) = 0.; // dr2/d(dudw1);
890  pm(1,2) = 0.; // dv2/d(dudw1);
891  pm(2,2) = dphi2ddudw1; // d(dudw2)/d(dudw1);
892  pm(3,2) = deta2ddudw1; // d(eta2)/d(dudw1);
893  pm(4,2) = 0.; // d(pinv2)/d(dudw1);
894 
895  pm(0,3) = 0.; // dr2/d(dvdw1);
896  pm(1,3) = 0.; // dv2/d(dvdw1);
897  pm(2,3) = dphi2ddvdw1; // d(phi2)/d(dvdw1);
898  pm(3,3) = deta2ddvdw1; // d(eta2)/d(dvdw1);
899  pm(4,3) = 0.; // d(pinv2)/d(dvdw1);
900 
901  pm(0,4) = 0.; // dr2/d(pinv1);
902  pm(1,4) = 0.; // dv2/d(pinv1);
903  pm(2,4) = 0.; // d(phi2)/d(pinv1);
904  pm(3,4) = 0.; // d(eta2)/d(pinv1);
905  pm(4,4) = 1.; // d(pinv2)/d(pinv1);
906  }
907 
908  // Update track vector.
909 
910  vec(0) = 0.;
911  vec(1) = 0.;
912  vec(2) = phid2;
913  vec(3) = eta2;
914 
915  // Done (success).
916 
917  return true;
918  }
TDirectory * dir
Definition: macro.C:5
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
bool trkf::PropYZLine::transformYZLine ( double  phi1,
double  phi2,
TrackVector vec,
Surface::TrackDirection dir,
TrackMatrix prop_matrix 
) const
private

Transform yz line -> yz line.

The following methods transform the track parameters from initial surface to SurfYZLine origin surface, and generate a propagation matrix. The first group of function parameters are the orientation surface parameters of the initial surface. The second group of function parameters are the orientation parameters of the of the destination surface. The origin parameters of the destination surface are assumed to match the position of the track.

Definition at line 376 of file PropYZLine.cxx.

Referenced by clone(), and origin_vec_prop().

380  {
381  // Calculate surface transcendental functions.
382 
383  double sindphi = std::sin(phi2 - phi1);
384  double cosdphi = std::cos(phi2 - phi1);
385 
386  // Get the initial track parameters.
387 
388  double r1 = vec(0);
389  double phid1 = vec(2);
390  double eta1 = vec(3);
391 
392  // Calculate elements of rotation matrix from initial coordinate
393  // system to destination coordinte system.
394 
395  double rvv = cosdphi;
396  double rvw = sindphi;
397 
398  double rwv = -sindphi;
399  double rww = cosdphi;
400 
401  // Calculate track transcendental functions.
402 
403  double sinphid1 = std::sin(phid1);
404  double cosphid1 = std::cos(phid1);
405  double sh1 = 1. / std::cosh(eta1); // sech(eta1)
406  double th1 = std::tanh(eta1);
407 
408  // Calculate initial position in Cartesian coordinates.
409 
410  double u1 = -r1 * sinphid1;
411  double w1 = r1 * cosphid1;
412 
413  // Calculate direction in destination coordinate system.
414 
415  double du2 = sh1*cosphid1;
416  double dv2 = th1*cosdphi + sh1*sinphid1*sindphi;
417  double dw2 = -th1*sindphi + sh1*sinphid1*cosdphi;
418  double duw2 = std::hypot(du2, dw2);
419 
420  // Calculate final direction track parameters.
421 
422  double phid2 = atan2(dw2, du2);
423  double eta2 = std::asinh(dv2 / duw2);
424 
425  // Update propagation matrix (if requested).
426 
427  if(prop_matrix != 0) {
428  TrackMatrix& pm = *prop_matrix;
429  pm.resize(vec.size(), vec.size(), false);
430 
431  // Calculate partial derivatives.
432 
433  // Partials of initial positions and directions wrt initial t.p.'s.
434 
435  double du1dr1 = -sinphid1;
436  double du1dphi1 = -w1;
437 
438  double dw1dr1 = cosphid1;
439  double dw1dphi1 = u1;
440 
441  double ddu1dphi1 = -sinphid1*sh1;
442  double ddu1deta1 = -cosphid1*sh1*th1;
443 
444  double ddv1deta1 = sh1*sh1;
445 
446  double ddw1dphi1 = cosphid1*sh1;
447  double ddw1deta1 = -sinphid1*sh1*th1;
448 
449  // Rotate partials to destination coordinate system.
450 
451  double du2dr1 = du1dr1;
452  double dv2dr1 = rvw*dw1dr1;
453  double dw2dr1 = rww*dw1dr1;
454 
455  double dv2dv1 = rvv;
456  double dw2dv1 = rwv;
457 
458  double du2dphi1 = du1dphi1;
459  double dv2dphi1 = rvw*dw1dphi1;
460  double dw2dphi1 = rww*dw1dphi1;
461 
462  double ddu2dphi1 = ddu1dphi1;
463  double ddv2dphi1 = rvw*ddw1dphi1;
464  double ddw2dphi1 = rww*ddw1dphi1;
465 
466  double ddu2deta1 = ddu1deta1;
467  double ddv2deta1 = rvv*ddv1deta1 + rvw*ddw1deta1;
468  double ddw2deta1 = rwv*ddv1deta1 + rww*ddw1deta1;
469 
470  // Partials of final t.p. wrt final position and direction.
471 
472  double dr2du2 = -dw2/duw2;
473  double dr2dw2 = du2/duw2;
474 
475  double dphi2ddu2 = -dw2/(duw2*duw2);
476  double dphi2ddw2 = du2/(duw2*duw2);
477 
478  double deta2ddv2 = 1./(duw2*duw2);
479 
480  // Partials of final t.p. wrt initial t.p.
481 
482  double dr2dr1 = dr2du2*du2dr1 + dr2dw2*dw2dr1;
483  double dr2dv1 = dr2dw2*dw2dv1;
484  double dr2dphi1 = dr2du2*du2dphi1 + dr2dw2*dw2dphi1;
485 
486  double dphi2dphi1 = dphi2ddu2*ddu2dphi1 + dphi2ddw2*ddw2dphi1;
487  double dphi2deta1 = dphi2ddu2*ddu2deta1 + dphi2ddw2*ddw2deta1;
488 
489  double deta2dphi1 = deta2ddv2*ddv2dphi1;
490  double deta2deta1 = deta2ddv2*ddv2deta1;
491 
492  // We still need to calculate the corretion due to the dependence of the
493  // propagation distance on the initial track parameters. This correction is
494  // needed even though the actual propagation distance is zero.
495 
496  // This correction only effects the v track parameter, since the v parameter
497  // the only parameter that actually dependes on the propagation distance.
498 
499  // Partials of propagation distance wrt position and direction in the destination
500  // coordinate system.
501 
502  double dsdu2 = -du2/(duw2*duw2);
503  double dsdw2 = -dw2/(duw2*duw2);
504 
505  // Partials of propagation distance wrt initial t.p.
506 
507  double dsdr1 = dsdu2*du2dr1 + dsdw2*dw2dr1;
508  double dsdv1 = dsdw2*dw2dv1;
509  double dsdphi1 = dsdu2*du2dphi1 + dsdw2*dw2dphi1;
510 
511  // Calculate correction to v parameter partials wrt initial t.p. due to path length.
512 
513  dv2dr1 += dv2*dsdr1;
514  dv2dv1 += dv2*dsdv1;
515  dv2dphi1 += dv2*dsdphi1;
516 
517  // Fill derivative matrix.
518 
519  pm(0,0) = dr2dr1; // dr2/dr1
520  pm(1,0) = dv2dr1; // dv2/dr1
521  pm(2,0) = 0.; // d(phi2)/dr1
522  pm(3,0) = 0.; // d(eta2)/dr1
523  pm(4,0) = 0.; // d(pinv2)/dr1
524 
525  pm(0,1) = dr2dv1; // dr2/dv1
526  pm(1,1) = dv2dv1; // dv2/dv1
527  pm(2,1) = 0.; // d(phi2)/dv1
528  pm(3,1) = 0.; // d(eta2)/dv1
529  pm(4,1) = 0.; // d(pinv2)/dv1
530 
531  pm(0,2) = dr2dphi1; // dr2/d(phi1);
532  pm(1,2) = dv2dphi1; // dv2/d(phi1);
533  pm(2,2) = dphi2dphi1; // d(phi2)/d(phi1);
534  pm(3,2) = deta2dphi1; // d(eta2)/d(phi1);
535  pm(4,2) = 0.; // d(pinv2)/d(phi1);
536 
537  pm(0,3) = 0.; // dr2/d(eta1);
538  pm(1,3) = 0.; // dv2/d(eta1);
539  pm(2,3) = dphi2deta1; // d(phi2)/d(eta1);
540  pm(3,3) = deta2deta1; // d(eta2)/d(eta1);
541  pm(4,3) = 0.; // d(pinv2)/d(eta1);
542 
543  pm(0,4) = 0.; // dr2/d(pinv1);
544  pm(1,4) = 0.; // dv2/d(pinv1);
545  pm(2,4) = 0.; // d(phi2)/d(pinv1);
546  pm(3,4) = 0.; // d(eta2)/d(pinv1);
547  pm(4,4) = 1.; // d(pinv2)/d(pinv1);
548  }
549 
550  // Update track vector.
551 
552  vec(0) = 0.;
553  vec(1) = 0.;
554  vec(2) = phid2;
555  vec(3) = eta2;
556 
557  // Done (success).
558 
559  return true;
560  }
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
bool trkf::PropYZLine::transformYZPlane ( double  phi1,
double  phi2,
TrackVector vec,
Surface::TrackDirection dir,
TrackMatrix prop_matrix 
) const
private

Transform yz plane -> yz line.

Definition at line 564 of file PropYZLine.cxx.

References trkf::Surface::BACKWARD, and trkf::Surface::FORWARD.

Referenced by clone(), and origin_vec_prop().

568  {
569  // Calculate surface transcendental functions.
570 
571  double sindphi = std::sin(phi2 - phi1);
572  double cosdphi = std::cos(phi2 - phi1);
573 
574  // Get the initial track parameters.
575 
576  double dudw1 = vec(2);
577  double dvdw1 = vec(3);
578 
579  // Make sure initial track has a valid direction.
580 
581  double dirf = 1.;
582  if(dir == Surface::BACKWARD)
583  dirf = -1.;
584  else if(dir != Surface::FORWARD)
585  return false;
586 
587  // Calculate elements of rotation matrix from initial coordinate
588  // system to destination coordinte system.
589 
590  double rvv = cosdphi;
591  double rvw = sindphi;
592 
593  double rwv = -sindphi;
594  double rww = cosdphi;
595 
596  // Calculate direction in the starting coordinate system.
597 
598  double dw1 = dirf / std::sqrt(1. + dudw1*dudw1 + dvdw1*dvdw1);
599  double du1 = dudw1 * dw1;
600  double dv1 = dvdw1 * dw1;
601 
602  // Rotate direction vector into destination coordinate system.
603 
604  double du2 = du1;
605  double dv2 = rvv*dv1 + rvw*dw1;
606  double dw2 = rwv*dv1 + rww*dw1;
607  double duw2 = std::hypot(du2, dw2);
608 
609  // Calculate final direction track parameters.
610 
611  double phid2 = atan2(dw2, du2);
612  double eta2 = std::asinh(dv2 / duw2);
613 
614  // Update propagation matrix (if requested).
615 
616  if(prop_matrix != 0) {
617  TrackMatrix& pm = *prop_matrix;
618  pm.resize(vec.size(), vec.size(), false);
619 
620  // Calculate partial derivatives.
621 
622  // Partials of initial positions and directions wrt initial t.p.'s.
623 
624  double ddu1ddudw1 = (1. + dvdw1*dvdw1) * dw1*dw1*dw1;
625  double ddu1ddvdw1 = -dudw1 * dvdw1 * dw1*dw1*dw1;
626 
627  double ddv1ddudw1 = -dudw1 * dvdw1 * dw1*dw1*dw1;
628  double ddv1ddvdw1 = (1. + dudw1*dudw1) * dw1*dw1*dw1;
629 
630  double ddw1ddudw1 = -dudw1 * dw1*dw1*dw1;
631  double ddw1ddvdw1 = -dvdw1 * dw1*dw1*dw1;
632 
633  // Rotate partials to destination coordinate system.
634 
635  double dv2dv1 = rvv;
636  double dw2dv1 = rwv;
637 
638  double ddu2ddudw1 = ddu1ddudw1;
639  double ddv2ddudw1 = rvv*ddv1ddudw1 + rvw*ddw1ddudw1;
640  double ddw2ddudw1 = rwv*ddv1ddudw1 + rww*ddw1ddudw1;
641 
642  double ddu2ddvdw1 = ddu1ddvdw1;
643  double ddv2ddvdw1 = rvv*ddv1ddvdw1 + rvw*ddw1ddvdw1;
644  double ddw2ddvdw1 = rwv*ddv1ddvdw1 + rww*ddw1ddvdw1;
645 
646  // Partials of final t.p. wrt final position and direction.
647 
648  double dr2du2 = -dw2/duw2;
649  double dr2dw2 = du2/duw2;
650 
651  double dphi2ddu2 = -dw2/(duw2*duw2);
652  double dphi2ddw2 = du2/(duw2*duw2);
653 
654  double deta2ddv2 = 1./(duw2*duw2);
655 
656  // Partials of final t.p. wrt initial t.p.
657 
658  double dr2du1 = dr2du2;
659  double dr2dv1 = dr2dw2*dw2dv1;
660 
661  double dphi2ddudw1 = dphi2ddu2*ddu2ddudw1 + dphi2ddw2*ddw2ddudw1;
662  double dphi2ddvdw1 = dphi2ddu2*ddu2ddvdw1 + dphi2ddw2*ddw2ddvdw1;
663 
664  double deta2ddudw1 = deta2ddv2*ddv2ddudw1;
665  double deta2ddvdw1 = deta2ddv2*ddv2ddvdw1;
666 
667  // We still need to calculate the corretion due to the dependence of the
668  // propagation distance on the initial track parameters. This correction is
669  // needed even though the actual propagation distance is zero.
670 
671  // This correction only effects the v track parameter, since the v parameter
672  // the only parameter that actually dependes on the propagation distance.
673 
674  // Partials of propagation distance wrt position and direction in the destination
675  // coordinate system.
676 
677  double dsdu2 = -du2/(duw2*duw2);
678  double dsdw2 = -dw2/(duw2*duw2);
679 
680  // Partials of propagation distance wrt initial t.p.
681 
682  double dsdu1 = dsdu2;
683  double dsdv1 = dsdw2*dw2dv1;
684 
685  // Calculate correction to v parameter partials wrt initial t.p. due to path length.
686 
687  double dv2du1 = dv2*dsdu1;
688  dv2dv1 += dv2*dsdv1;
689 
690  // Fill matrix.
691 
692  pm(0,0) = dr2du1; // dr2/du1
693  pm(1,0) = dv2du1; // dv2/du1
694  pm(2,0) = 0.; // d(phi2)/du1
695  pm(3,0) = 0.; // d(eta2)/du1
696  pm(4,0) = 0.; // d(pinv2)/du1
697 
698  pm(0,1) = dr2dv1; // dr2/dv1
699  pm(1,1) = dv2dv1; // dv2/dv1
700  pm(2,1) = 0.; // d(phi2)/dv1
701  pm(3,1) = 0.; // d(eta2)/dv1
702  pm(4,1) = 0.; // d(pinv2)/dv1
703 
704  pm(0,2) = 0.; // dr2/d(dudw1);
705  pm(1,2) = 0.; // dv2/d(dudw1);
706  pm(2,2) = dphi2ddudw1; // d(dudw2)/d(dudw1);
707  pm(3,2) = deta2ddudw1; // d(eta2)/d(dudw1);
708  pm(4,2) = 0.; // d(pinv2)/d(dudw1);
709 
710  pm(0,3) = 0.; // dr2/d(dvdw1);
711  pm(1,3) = 0.; // dv2/d(dvdw1);
712  pm(2,3) = dphi2ddvdw1; // d(phi2)/d(dvdw1);
713  pm(3,3) = deta2ddvdw1; // d(eta2)/d(dvdw1);
714  pm(4,3) = 0.; // d(pinv2)/d(dvdw1);
715 
716  pm(0,4) = 0.; // dr2/d(pinv1);
717  pm(1,4) = 0.; // dv2/d(pinv1);
718  pm(2,4) = 0.; // d(phi2)/d(pinv1);
719  pm(3,4) = 0.; // d(eta2)/d(pinv1);
720  pm(4,4) = 1.; // d(pinv2)/d(pinv1);
721  }
722 
723  // Update track vector.
724 
725  vec(0) = 0.;
726  vec(1) = 0.;
727  vec(2) = phid2;
728  vec(3) = eta2;
729 
730  // Done (success).
731 
732  return true;
733  }
TDirectory * dir
Definition: macro.C:5
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
boost::optional< double > trkf::Propagator::vec_prop ( KTrack trk,
const std::shared_ptr< const Surface > &  psurf,
PropDirection  dir,
bool  doDedx,
TrackMatrix prop_matrix = 0,
TrackError noise_matrix = 0 
) const
inherited

Propagate without error (long distance).

Propagate without error (long distance).

Arguments:

trk - Track to propagate. psurf - Destination surface. dir - Propagation direction (FORWARD, BACKWARD, or UNKNOWN). doDedx - dE/dx enable/disable flag. prop_matrix - Return propagation matrix if not null. noise_matrix - Return noise matrix if not null.

Returned value: Propagation distance & success flag.

This method calls virtual method short_vec_prop in steps of some maximum size.

Definition at line 52 of file Propagator.cxx.

References e, trkf::Propagator::fTcut, trkf::Propagator::getDoDedx(), trkf::KTrack::getMomentum(), trkf::KTrack::getPosition(), trkf::KTrack::getVector(), trkf::KTrack::Mass(), s, and trkf::Propagator::short_vec_prop().

Referenced by trkf::KalmanFilterAlg::fitMomentumMS(), trkf::Propagator::getInteractor(), trkf::Propagator::lin_prop(), trkf::KHitContainer::sort(), and trkf::KalmanFilterAlg::updateMomentum().

58  {
59  // Default result.
60 
61  auto result = boost::make_optional<double>(false, 0.);
62 
63  // Get the inverse momentum (assumed to be track parameter four).
64 
65  double pinv = trk.getVector()(4);
66 
67  // If dE/dx is not requested, or if inverse momentum is zero, then
68  // it is safe to propagate in one step. In this case, just pass
69  // the call to short_vec_prop with unlimited distance.
70 
71  bool dedx = getDoDedx() && doDedx;
72  if(!dedx || pinv == 0.)
73  result = short_vec_prop(trk, psurf, dir, dedx, prop_matrix, noise_matrix);
74 
75  else {
76 
77  // dE/dx is requested. In this case we limit the maximum
78  // propagation distance such that the kinetic energy of the
79  // particle should not change by more thatn 10%.
80 
81  // Get LAr service.
82 
83  auto const * detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
84 
85  // Initialize propagation matrix to unit matrix (if specified).
86 
87  int nvec = trk.getVector().size();
88  if(prop_matrix)
89  *prop_matrix = ublas::identity_matrix<TrackVector::value_type>(nvec);
90 
91  // Initialize noise matrix to zero matrix (if specified).
92 
93  if(noise_matrix) {
94  noise_matrix->resize(nvec, nvec, false);
95  noise_matrix->clear();
96  }
97 
98  // Remember the starting track.
99 
100  KTrack trk0(trk);
101 
102  // Make pointer variables pointing to local versions of the
103  // propagation and noise matrices, or null if not specified.
104 
105  TrackMatrix local_prop_matrix;
106  TrackMatrix* plocal_prop_matrix = (prop_matrix==0 ? 0 : &local_prop_matrix);
107  TrackError local_noise_matrix;
108  TrackError* plocal_noise_matrix = (noise_matrix==0 ? 0 : &local_noise_matrix);
109 
110  // Cumulative propagation distance.
111 
112  double s = 0.;
113 
114  // Begin stepping loop.
115  // We put a maximum iteration count to prevent infinite loops caused by
116  // floating point pathologies. The iteration count is large enough to reach
117  // any point in the tpc using the minimum step size (for a reasonable tpc).
118 
119  bool done = false;
120  int nitmax = 10000; // Maximum number of iterations.
121  int nit = 0; // Iteration count.
122  while(!done) {
123 
124  // If the iteration count exceeds the maximum, return failure.
125 
126  ++nit;
127  if(nit > nitmax) {
128  trk = trk0;
129  result = boost::optional<double>(false, 0.);
130  return result;
131  }
132 
133  // Estimate maximum step distance according to the above
134  // stated principle.
135 
136  pinv = trk.getVector()(4);
137  double mass = trk.Mass();
138  double p = 1./std::abs(pinv);
139  double e = std::hypot(p, mass);
140  double t = p*p / (e + mass);
141  double dedx = 0.001 * detprop->Eloss(p, mass, fTcut);
142  double smax = 0.1 * t / dedx;
143  if (smax <= 0.)
144  throw cet::exception("Propagator") << __func__ << ": maximum step " << smax << "\n";
145 
146  // Always allow a step of at least 0.3 cm (about one wire spacing).
147 
148  if(smax < 0.3)
149  smax = 0.3;
150 
151  // First do a test propagation (without dE/dx and errors) to
152  // find the distance to the destination surface.
153 
154  KTrack trktest(trk);
155  boost::optional<double> dist = short_vec_prop(trktest, psurf, dir, false, 0, 0);
156 
157  // If the test propagation failed, return failure.
158 
159  if(!dist) {
160  trk = trk0;
161  return dist;
162  }
163 
164  // Generate destionation surface for this step (either final
165  // destination, or some intermediate surface).
166 
167  std::shared_ptr<const Surface> pstep;
168  if(std::abs(*dist) <= smax) {
169  done = true;
170  pstep = psurf;
171  }
172  else {
173 
174  // Generate intermediate surface.
175  // First get point where track will intersect intermediate surface.
176 
177  double xyz0[3]; // Starting point.
178  trk.getPosition(xyz0);
179  double xyz1[3]; // Destination point.
180  trktest.getPosition(xyz1);
181  double frac = smax / std::abs(*dist);
182  double xyz[3]; // Intermediate point.
183  xyz[0] = xyz0[0] + frac * (xyz1[0] - xyz0[0]);
184  xyz[1] = xyz0[1] + frac * (xyz1[1] - xyz0[1]);
185  xyz[2] = xyz0[2] + frac * (xyz1[2] - xyz0[2]);
186 
187  // Choose orientation of intermediate surface perpendicular
188  // to track.
189 
190  double mom[3];
191  trk.getMomentum(mom);
192 
193  // Make intermediate surface object.
194 
195  pstep = std::shared_ptr<const Surface>(new SurfXYZPlane(xyz[0], xyz[1], xyz[2],
196  mom[0], mom[1], mom[2]));
197  }
198 
199  // Do the actual step propagation.
200 
201  dist = short_vec_prop(trk, pstep, dir, doDedx,
202  plocal_prop_matrix, plocal_noise_matrix);
203 
204  // If the step propagation failed, return failure.
205 
206  if(!dist) {
207  trk = trk0;
208  return dist;
209  }
210 
211  // Update cumulative propagation distance.
212 
213  s += *dist;
214 
215  // Update cumulative propagation matrix (left-multiply).
216 
217  if(prop_matrix != 0) {
218  TrackMatrix temp = prod(*plocal_prop_matrix, *prop_matrix);
219  *prop_matrix = temp;
220  }
221 
222  // Update cumulative noise matrix.
223 
224  if(noise_matrix != 0) {
225  TrackMatrix temp = prod(*noise_matrix, trans(*plocal_prop_matrix));
226  TrackMatrix temp2 = prod(*plocal_prop_matrix, temp);
227  *noise_matrix = ublas::symmetric_adaptor<TrackMatrix>(temp2);
228  *noise_matrix += *plocal_noise_matrix;
229  }
230  }
231 
232  // Set the final result (distance + success).
233 
234  result = boost::optional<double>(true, s);
235  }
236 
237  // Done.
238 
239  return result;
240  }
Float_t s
Definition: plot.C:23
KSymMatrix< 5 >::type TrackError
Track error matrix, dimension 5x5.
double fTcut
Maximum delta ray energy for dE/dx.
Definition: Propagator.h:163
bool getDoDedx() const
Definition: Propagator.h:103
TDirectory * dir
Definition: macro.C:5
KMatrix< 5, 5 >::type TrackMatrix
General 5x5 matrix.
Float_t e
Definition: plot.C:34
virtual boost::optional< double > short_vec_prop(KTrack &trk, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const =0
Propagate without error (short distance).
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33

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