CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Friends
CachedTrajectory Class Reference

#include <CachedTrajectory.h>

Public Member Functions

const std::vector< SteppingHelixStateInfo > & getEcalTrajectory () const
 
const std::vector< SteppingHelixStateInfo > & getHcalTrajectory () const
 
const std::vector< SteppingHelixStateInfo > & getHOTrajectory () const
 
const std::vector< SteppingHelixStateInfo > & getPreshowerTrajectory () const
 

Protected Member Functions

std::pair< float, float > delta (const double &theta1, const double &theta2, const double &phi1, const double &phi2)
 
float distance (const Plane *plane, int index)
 

Static Protected Member Functions

static int sign (float number)
 

Protected Attributes

std::vector< SteppingHelixStateInfoecalTrajectory_
 
std::deque< SteppingHelixStateInfofullTrajectory_
 
bool fullTrajectoryFilled_
 
std::vector< SteppingHelixStateInfohcalTrajectory_
 
float HOmaxRho_
 
float HOmaxZ_
 
std::vector< SteppingHelixStateInfohoTrajectory_
 
float maxRho_
 
float maxZ_
 
float minRho_
 
float minZ_
 
std::vector< SteppingHelixStateInfopreshowerTrajectory_
 
const Propagatorpropagator_
 
SteppingHelixStateInfo stateAtIP_
 
float step_
 
std::vector< GlobalPointwideEcalTrajectory_
 
std::vector< GlobalPointwideHcalTrajectory_
 
std::vector< GlobalPointwideHOTrajectory_
 

Private Types

enum  TrajectorType { IpToEcal, IpToHcal, IpToHO, FullTrajectory }
 
enum  WideTrajectoryType { Ecal, Hcal, HO }
 

Private Member Functions

 CachedTrajectory ()
 
void findEcalTrajectory (const FiducialVolume &)
 
void findHcalTrajectory (const FiducialVolume &)
 
void findHOTrajectory (const FiducialVolume &)
 
void findPreshowerTrajectory (const FiducialVolume &)
 
SteppingHelixStateInfo getInnerState ()
 
SteppingHelixStateInfo getOuterState ()
 
float getPropagationStep () const
 
SteppingHelixStateInfo getStateAtEcal ()
 
SteppingHelixStateInfo getStateAtHcal ()
 
SteppingHelixStateInfo getStateAtHO ()
 
SteppingHelixStateInfo getStateAtPreshower ()
 
void getTrajectory (std::vector< SteppingHelixStateInfo > &, const FiducialVolume &, int steps=4)
 
std::vector< GlobalPoint > * getWideTrajectory (const std::vector< SteppingHelixStateInfo > &, WideTrajectoryType)
 
TrajectoryStateOnSurface propagate (const Plane *plane)
 get fast to a given DetId surface using cached trajectory More...
 
void propagate (SteppingHelixStateInfo &state, const Cylinder &cylinder)
 
void propagate (SteppingHelixStateInfo &state, const Plane &plane)
 
bool propagateAll (const SteppingHelixStateInfo &initialState)
 propagate through the whole detector, returns true if successful More...
 
void propagateForward (SteppingHelixStateInfo &state, float distance)
 
void reset_trajectory ()
 
void setMaxDetectorLength (float l=2200.)
 
void setMaxDetectorRadius (float r=800.)
 
void setMaxHOLength (float l=2200.)
 
void setMaxHORadius (float r=800.)
 
void setMinDetectorLength (float l=0.)
 
void setMinDetectorRadius (float r=0.)
 
void setPropagationStep (float s=20.)
 
void setPropagator (const Propagator *ptr)
 
void setStateAtIP (const SteppingHelixStateInfo &state)
 
std::pair< float, float > trajectoryDelta (TrajectorType)
 

Friends

std::vector< SteppingHelixStateInfopropagateThoughFromIP (const SteppingHelixStateInfo &state, const Propagator *ptr, const FiducialVolume &volume, int nsteps, float step, float minR, float minZ, float maxR, float maxZ)
 
class TrackDetectorAssociator
 

Detailed Description

Definition at line 44 of file CachedTrajectory.h.

Member Enumeration Documentation

◆ TrajectorType

Enumerator
IpToEcal 
IpToHcal 
IpToHO 
FullTrajectory 

Definition at line 64 of file CachedTrajectory.h.

◆ WideTrajectoryType

Enumerator
Ecal 
Hcal 
HO 

Definition at line 65 of file CachedTrajectory.h.

65 { Ecal, Hcal, HO };

Constructor & Destructor Documentation

◆ CachedTrajectory()

CachedTrajectory::CachedTrajectory ( )
private

Member Function Documentation

◆ delta()

std::pair< float, float > CachedTrajectory::delta ( const double &  theta1,
const double &  theta2,
const double &  phi1,
const double &  phi2 
)
protected

Definition at line 319 of file CachedTrajectory.cc.

322  {
323  std::pair<float, float> result(theta2 - theta1, phi2 - phi1);
324  // this won't work for loopers, since deltaPhi cannot be larger than Pi.
325  if (fabs(result.second) > 2 * M_PI - fabs(result.second)) {
326  if (result.second > 0)
327  result.second -= 2 * M_PI;
328  else
329  result.second += 2 * M_PI;
330  }
331  return result;
332 }

References M_PI, and mps_fire::result.

Referenced by trajectoryDelta().

◆ distance()

float CachedTrajectory::distance ( const Plane plane,
int  index 
)
inlineprotected

Definition at line 136 of file CachedTrajectory.h.

136  {
137  if (index < 0 || fullTrajectory_.empty() || (unsigned int)index >= fullTrajectory_.size())
138  return 0;
139  return plane->localZ(fullTrajectory_[index].position());
140  }

References fullTrajectory_, createfilelist::int, and position.

Referenced by propagate(), and propagateForward().

◆ findEcalTrajectory()

void CachedTrajectory::findEcalTrajectory ( const FiducialVolume volume)
private

Definition at line 443 of file CachedTrajectory.cc.

443  {
444  LogTrace("TrackAssociator") << "getting trajectory in ECAL";
445  getTrajectory(ecalTrajectory_, volume, 4);
446  LogTrace("TrackAssociator") << "# of points in ECAL trajectory:" << ecalTrajectory_.size();
447 }

References ecalTrajectory_, getTrajectory(), and LogTrace.

◆ findHcalTrajectory()

void CachedTrajectory::findHcalTrajectory ( const FiducialVolume volume)
private

Definition at line 461 of file CachedTrajectory.cc.

461  {
462  LogTrace("TrackAssociator") << "getting trajectory in HCAL";
463  getTrajectory(hcalTrajectory_, volume, 4); // more steps to account for different depth
464  LogTrace("TrackAssociator") << "# of points in HCAL trajectory:" << hcalTrajectory_.size();
465 }

References getTrajectory(), hcalTrajectory_, and LogTrace.

◆ findHOTrajectory()

void CachedTrajectory::findHOTrajectory ( const FiducialVolume volume)
private

Definition at line 469 of file CachedTrajectory.cc.

469  {
470  LogTrace("TrackAssociator") << "getting trajectory in HO";
471  getTrajectory(hoTrajectory_, volume, 2);
472  LogTrace("TrackAssociator") << "# of points in HO trajectory:" << hoTrajectory_.size();
473 }

References getTrajectory(), hoTrajectory_, and LogTrace.

◆ findPreshowerTrajectory()

void CachedTrajectory::findPreshowerTrajectory ( const FiducialVolume volume)
private

Definition at line 449 of file CachedTrajectory.cc.

449  {
450  LogTrace("TrackAssociator") << "getting trajectory in Preshower";
452  LogTrace("TrackAssociator") << "# of points in Preshower trajectory:" << preshowerTrajectory_.size();
453 }

References getTrajectory(), LogTrace, and preshowerTrajectory_.

◆ getEcalTrajectory()

const std::vector< SteppingHelixStateInfo > & CachedTrajectory::getEcalTrajectory ( ) const

Definition at line 455 of file CachedTrajectory.cc.

455 { return ecalTrajectory_; }

References ecalTrajectory_.

◆ getHcalTrajectory()

const std::vector< SteppingHelixStateInfo > & CachedTrajectory::getHcalTrajectory ( ) const

Definition at line 467 of file CachedTrajectory.cc.

467 { return hcalTrajectory_; }

References hcalTrajectory_.

◆ getHOTrajectory()

const std::vector< SteppingHelixStateInfo > & CachedTrajectory::getHOTrajectory ( ) const

Definition at line 475 of file CachedTrajectory.cc.

475 { return hoTrajectory_; }

References hoTrajectory_.

◆ getInnerState()

SteppingHelixStateInfo CachedTrajectory::getInnerState ( )
private

Definition at line 598 of file CachedTrajectory.cc.

598  {
599  if (fullTrajectory_.empty())
600  return SteppingHelixStateInfo();
601  else
602  return fullTrajectory_.front();
603 }

References fullTrajectory_.

◆ getOuterState()

SteppingHelixStateInfo CachedTrajectory::getOuterState ( )
private

Definition at line 605 of file CachedTrajectory.cc.

605  {
606  if (fullTrajectory_.empty())
607  return SteppingHelixStateInfo();
608  else
609  return fullTrajectory_.back();
610 }

References fullTrajectory_.

◆ getPreshowerTrajectory()

const std::vector< SteppingHelixStateInfo > & CachedTrajectory::getPreshowerTrajectory ( ) const

Definition at line 457 of file CachedTrajectory.cc.

457  {
458  return preshowerTrajectory_;
459 }

References preshowerTrajectory_.

◆ getPropagationStep()

float CachedTrajectory::getPropagationStep ( ) const
inlineprivate

Definition at line 121 of file CachedTrajectory.h.

121 { return step_; }

References step_.

◆ getStateAtEcal()

SteppingHelixStateInfo CachedTrajectory::getStateAtEcal ( )
private

Definition at line 570 of file CachedTrajectory.cc.

570  {
571  if (ecalTrajectory_.empty())
572  return SteppingHelixStateInfo();
573  else
574  return ecalTrajectory_.front();
575 }

References ecalTrajectory_.

◆ getStateAtHcal()

SteppingHelixStateInfo CachedTrajectory::getStateAtHcal ( )
private

Definition at line 584 of file CachedTrajectory.cc.

584  {
585  if (hcalTrajectory_.empty())
586  return SteppingHelixStateInfo();
587  else
588  return hcalTrajectory_.front();
589 }

References hcalTrajectory_.

◆ getStateAtHO()

SteppingHelixStateInfo CachedTrajectory::getStateAtHO ( )
private

Definition at line 591 of file CachedTrajectory.cc.

591  {
592  if (hoTrajectory_.empty())
593  return SteppingHelixStateInfo();
594  else
595  return hoTrajectory_.front();
596 }

References hoTrajectory_.

◆ getStateAtPreshower()

SteppingHelixStateInfo CachedTrajectory::getStateAtPreshower ( )
private

Definition at line 577 of file CachedTrajectory.cc.

577  {
578  if (preshowerTrajectory_.empty())
579  return SteppingHelixStateInfo();
580  else
581  return preshowerTrajectory_.front();
582 }

References preshowerTrajectory_.

◆ getTrajectory()

void CachedTrajectory::getTrajectory ( std::vector< SteppingHelixStateInfo > &  trajectory,
const FiducialVolume volume,
int  steps = 4 
)
private

get a set of points representing the trajectory between two cylinders of radius R1 and R2 and length L1 and L2. Parameter steps defines maximal number of steps in the detector.

Definition at line 334 of file CachedTrajectory.cc.

336  {
338  throw cms::Exception("FatalError") << "trajectory is not defined yet. Please use propagateAll first.";
339  if (fullTrajectory_.empty()) {
340  LogTrace("TrackAssociator") << "Trajectory is empty. Move on";
341  return;
342  }
343  if (!volume.isValid()) {
344  LogTrace("TrackAssociator") << "no trajectory is expected to be found since the fiducial volume is not valid";
345  return;
346  }
347  double step = std::max(volume.maxR() - volume.minR(), volume.maxZ() - volume.minZ()) / steps;
348 
349  int closestPointOnLeft = -1;
350 
351  // check whether the trajectory crossed the region
352  if (!((fullTrajectory_.front().position().perp() < volume.maxR() &&
353  fabs(fullTrajectory_.front().position().z()) < volume.maxZ()) &&
354  (fullTrajectory_.back().position().perp() > volume.minR() ||
355  fabs(fullTrajectory_.back().position().z()) > volume.minZ()))) {
356  LogTrace("TrackAssociator") << "Track didn't cross the region (R1,R2,L1,L2): " << volume.minR() << ", "
357  << volume.maxR() << ", " << volume.minZ() << ", " << volume.maxZ();
358  return;
359  }
360 
361  // get distance along momentum to the surface.
362 
363  // the following code can be made faster, but it'll hardly be a significant improvement
364  // simplifications:
365  // 1) direct loop over stored trajectory points instead of some sort
366  // of fast root search (Newton method)
367  // 2) propagate from the closest point outside the region with the
368  // requested step ignoring stored trajectory points.
369  double dZ(-1.);
370  double dR(-1.);
371  int firstPointInside(-1);
372  for (unsigned int i = 0; i < fullTrajectory_.size(); i++) {
373  // LogTrace("TrackAssociator") << "Trajectory info (i,perp,r1,r2,z,z1,z2): " << i << ", " << fullTrajectory_[i].position().perp() <<
374  // ", " << volume.minR() << ", " << volume.maxR() << ", " << fullTrajectory_[i].position().z() << ", " << volume.minZ() << ", " <<
375  // volume.maxZ() << ", " << closestPointOnLeft;
376  dR = fullTrajectory_[i].position().perp() - volume.minR();
377  dZ = fabs(fullTrajectory_[i].position().z()) - volume.minZ();
378  if (dR > 0 || dZ > 0) {
379  if (i > 0) {
380  firstPointInside = i;
381  closestPointOnLeft = i - 1;
382  } else {
383  firstPointInside = 0;
384  closestPointOnLeft = 0;
385  }
386  break;
387  }
388  }
389  if (closestPointOnLeft == -1)
390  throw cms::Exception("FatalError") << "This shouls never happen - internal logic error";
391 
392  SteppingHelixStateInfo currentState(fullTrajectory_[closestPointOnLeft]);
393  if (currentState.position().x() * currentState.momentum().x() +
394  currentState.position().y() * currentState.momentum().y() +
395  currentState.position().z() * currentState.momentum().z() <
396  0)
397  step = -step;
398 
399  // propagate to the inner surface of the active volume
400 
401  if (firstPointInside != closestPointOnLeft) {
402  if (dR > 0) {
405  propagate(currentState, *barrel);
406  } else {
408  Plane::build(Plane::PositionType(0, 0, currentState.position().z() > 0 ? volume.minZ() : -volume.minZ()),
410  propagate(currentState, *endcap);
411  }
412  if (currentState.isValid())
413  trajectory.push_back(currentState);
414  } else
415  LogTrace("TrackAssociator") << "Weird message\n";
416 
417  while (currentState.isValid() && currentState.position().perp() < volume.maxR() &&
418  fabs(currentState.position().z()) < volume.maxZ()) {
419  propagateForward(currentState, step);
420  if (!currentState.isValid()) {
421  LogTrace("TrackAssociator") << "Failed to propagate the track; moving on\n";
422  break;
423  }
424  // LogTrace("TrackAssociator") << "New state (perp, z): " << currentState.position().perp() << ", " << currentState.position().z();
425  //if ( ( currentState.position().perp() < volume.maxR() && fabs(currentState.position().z()) < volume.maxZ() ) &&
426  // ( currentState.position().perp()-volume.minR() > 0 || fabs(currentState.position().z()) - volume.minZ() >0 ) )
427  trajectory.push_back(currentState);
428  }
429 }

References Reference_intrackfit_cff::barrel, Plane::build(), Cylinder::build(), HGC3DClusterGenMatchSelector_cfi::dR, makeMuonMisalignmentScenario::endcap, Exception, fullTrajectory_, fullTrajectoryFilled_, mps_fire::i, FiducialVolume::isValid(), SteppingHelixStateInfo::isValid(), LogTrace, SiStripPI::max, FiducialVolume::maxR(), FiducialVolume::maxZ(), FiducialVolume::minR(), FiducialVolume::minZ(), SteppingHelixStateInfo::momentum(), PV3DBase< T, PVType, FrameType >::perp(), SteppingHelixStateInfo::position(), position, propagate(), propagateForward(), customisers::steps, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by findEcalTrajectory(), findHcalTrajectory(), findHOTrajectory(), findPreshowerTrajectory(), and propagateThoughFromIP().

◆ getWideTrajectory()

std::vector< GlobalPoint > * CachedTrajectory::getWideTrajectory ( const std::vector< SteppingHelixStateInfo > &  states,
WideTrajectoryType  wideTrajectoryType 
)
private

Definition at line 477 of file CachedTrajectory.cc.

478  {
479  std::vector<GlobalPoint>* wideTrajectory = nullptr;
480  switch (wideTrajectoryType) {
481  case Ecal:
482  LogTrace("TrackAssociator") << "Filling ellipses in Ecal trajectory";
483  wideTrajectory = &wideEcalTrajectory_;
484  break;
485  case Hcal:
486  LogTrace("TrackAssociator") << "Filling ellipses in Hcal trajectory";
487  wideTrajectory = &wideHcalTrajectory_;
488  break;
489  case HO:
490  LogTrace("TrackAssociator") << "Filling ellipses in HO trajectory";
491  wideTrajectory = &wideHOTrajectory_;
492  break;
493  }
494  if (!wideTrajectory)
495  return nullptr;
496 
497  for (std::vector<SteppingHelixStateInfo>::const_iterator state = states.begin(); state != states.end(); state++) {
498  // defined a normal plane wrt the particle trajectory direction
499  // let's hope that I computed the rotation matrix correctly.
500  GlobalVector vector(state->momentum().unit());
501  float r21 = 0;
502  float r22 = vector.z() / sqrt(1 - pow(vector.x(), 2));
503  float r23 = -vector.y() / sqrt(1 - pow(vector.x(), 2));
504  float r31 = vector.x();
505  float r32 = vector.y();
506  float r33 = vector.z();
507  float r11 = r22 * r33 - r23 * r32;
508  float r12 = r23 * r31;
509  float r13 = -r22 * r31;
510 
511  Plane::RotationType rotation(r11, r12, r13, r21, r22, r23, r31, r32, r33);
512  Plane* target = new Plane(state->position(), rotation);
513 
514  TrajectoryStateOnSurface tsos = state->getStateOnSurface(*target);
515 
516  if (!tsos.isValid()) {
517  LogTrace("TrackAssociator") << "[getWideTrajectory] TSOS not valid";
518  continue;
519  }
520  if (!tsos.hasError()) {
521  LogTrace("TrackAssociator") << "[getWideTrajectory] TSOS does not have Errors";
522  continue;
523  }
524  LocalError localErr = tsos.localError().positionError();
525  localErr.scale(2); // get the 2 sigma ellipse
526  float xx = localErr.xx();
527  float xy = localErr.xy();
528  float yy = localErr.yy();
529 
530  float denom = yy - xx;
531  float phi = 0.;
532  if (xy == 0 && denom == 0)
533  phi = M_PI_4;
534  else
535  phi = 0.5 * atan2(2. * xy, denom); // angle of MAJOR axis
536  // Unrotate the error ellipse to get the semimajor and minor axes. Then place points on
537  // the endpoints of semiminor an seminajor axes on original(rotated) error ellipse.
538  LocalError rotErr = localErr.rotate(-phi); // xy covariance of rotErr should be zero
539  float semi1 = sqrt(rotErr.xx());
540  float semi2 = sqrt(rotErr.yy());
541 
542  // Just use one point if the ellipse is small
543  // if(semi1 < 0.1 && semi2 < 0.1) {
544  // LogTrace("TrackAssociator") << "[getWideTrajectory] Error ellipse is small, using one trajectory point";
545  // wideTrajectory->push_back(state->position());
546  // continue;
547  // }
548 
549  Local2DPoint bounds[4];
550  bounds[0] = Local2DPoint(semi1 * cos(phi), semi1 * sin(phi));
551  bounds[1] = Local2DPoint(semi1 * cos(phi + M_PI), semi1 * sin(phi + M_PI));
552  phi += M_PI_2; // add pi/2 for the semi2 axis
553  bounds[2] = Local2DPoint(semi2 * cos(phi), semi2 * sin(phi));
554  bounds[3] = Local2DPoint(semi2 * cos(phi + M_PI), semi2 * sin(phi + M_PI));
555 
556  // LogTrace("TrackAssociator") << "Axes " << semi1 <<","<< semi2 <<" phi "<< phi;
557  // LogTrace("TrackAssociator") << "Local error ellipse: " << bounds[0] << bounds[1] << bounds[2] << bounds[3];
558 
559  wideTrajectory->push_back(state->position());
560  for (int index = 0; index < 4; ++index)
561  wideTrajectory->push_back(target->toGlobal(bounds[index]));
562 
563  //LogTrace("TrackAssociator") <<"Global error ellipse: (" << target->toGlobal(bounds[0]) <<","<< target->toGlobal(bounds[1])
564  // <<","<< target->toGlobal(bounds[2]) <<","<< target->toGlobal(bounds[3]) <<","<<state->position() <<")";
565  }
566 
567  return wideTrajectory;
568 }

References funct::cos(), makePileupJSON::denom, Ecal, Hcal, HO, LogTrace, M_PI, M_PI_2, phi, funct::pow(), dttmaxenums::r32, LocalError::rotate(), idealTransformation::rotation, LocalError::scale(), funct::sin(), mathSSE::sqrt(), filterCSVwithJSON::target, trackerHitRTTI::vector, wideEcalTrajectory_, wideHcalTrajectory_, wideHOTrajectory_, geometryCSVtoXML::xx, LocalError::xx(), geometryCSVtoXML::xy, LocalError::xy(), geometryCSVtoXML::yy, and LocalError::yy().

◆ propagate() [1/3]

TrajectoryStateOnSurface CachedTrajectory::propagate ( const Plane plane)
private

get fast to a given DetId surface using cached trajectory

Definition at line 179 of file CachedTrajectory.cc.

179  {
180  // TimerStack timers(TimerStack::Disableable);
181  // timers.benchmark("CachedTrajectory::propagate::benchmark");
182  // timers.push("CachedTrajectory::propagate",TimerStack::FastMonitoring);
183  // timers.push("CachedTrajectory::propagate::findClosestPoint",TimerStack::FastMonitoring);
184 
185  // Assume that all points along the trajectory are equally spread out.
186  // For simplication assume that the trajectory is just a straight
187  // line and find a point closest to the target plane. Propagate to
188  // the plane from the point.
189 
190  const float matchingDistance = 1;
191  // find the closest point to the plane
192  int leftIndex = 0;
193  int rightIndex = fullTrajectory_.size() - 1;
194  int closestPointOnLeft = 0;
195 
196  // check whether the trajectory crossed the plane (signs should be different)
197  if (sign(distance(plane, leftIndex)) * sign(distance(plane, rightIndex)) != -1) {
198  LogTrace("TrackAssociator") << "Track didn't cross the plane:\n\tleft distance: " << distance(plane, leftIndex)
199  << "\n\tright distance: " << distance(plane, rightIndex);
200  return TrajectoryStateOnSurface();
201  }
202 
203  while (leftIndex + 1 < rightIndex) {
204  closestPointOnLeft = int((leftIndex + rightIndex) / 2);
205  float dist = distance(plane, closestPointOnLeft);
206  // LogTrace("TrackAssociator") << "Closest point on left: " << closestPointOnLeft << "\n"
207  // << "Distance to the plane: " << dist;
208  if (fabs(dist) < matchingDistance) {
209  // found close match, verify that we are on the left side
210  if (closestPointOnLeft > 0 && sign(distance(plane, closestPointOnLeft - 1)) * dist == -1)
211  closestPointOnLeft--;
212  break;
213  }
214 
215  // check where is the plane
216  if (sign(distance(plane, leftIndex) * dist) == -1)
217  rightIndex = closestPointOnLeft;
218  else
219  leftIndex = closestPointOnLeft;
220 
221  // LogTrace("TrackAssociator") << "Distance on left: " << distance(plane, leftIndex) << "\n"
222  // << "Distance to closest point: " << distance(plane, closestPointOnLeft) << "\n"
223  // << "Left index: " << leftIndex << "\n"
224  // << "Right index: " << rightIndex;
225  }
226  LogTrace("TrackAssociator") << "closestPointOnLeft: " << closestPointOnLeft << "\n\ttrajectory point (z,R,eta,phi): "
227  << fullTrajectory_[closestPointOnLeft].position().z() << ", "
228  << fullTrajectory_[closestPointOnLeft].position().perp() << " , "
229  << fullTrajectory_[closestPointOnLeft].position().eta() << " , "
230  << fullTrajectory_[closestPointOnLeft].position().phi()
231  << "\n\tplane center (z,R,eta,phi): " << plane->position().z() << ", "
232  << plane->position().perp() << " , " << plane->position().eta() << " , "
233  << plane->position().phi();
234 
235  // propagate to the plane
236  // timers.pop_and_push("CachedTrajectory::propagate::localPropagation",TimerStack::FastMonitoring);
237  if (const SteppingHelixPropagator* shp = dynamic_cast<const SteppingHelixPropagator*>(propagator_)) {
239  try {
240  shp->propagate(fullTrajectory_[closestPointOnLeft], *plane, state);
241  } catch (cms::Exception& ex) {
242  edm::LogWarning("TrackAssociator") << "Caught exception " << ex.category() << ": " << ex.explainSelf();
243  edm::LogWarning("TrackAssociator") << "An exception is caught during the track propagation\n"
244  << state.momentum().x() << ", " << state.momentum().y() << ", "
245  << state.momentum().z();
246  return TrajectoryStateOnSurface();
247  }
248  return state.getStateOnSurface(*plane);
249  } else {
251  fullTrajectory_[closestPointOnLeft].getFreeState(fts);
252  return propagator_->propagate(fts, *plane);
253  }
254 }

References cms::Exception::category(), distance(), PV3DBase< T, PVType, FrameType >::eta(), cms::Exception::explainSelf(), fullTrajectory_, createfilelist::int, LogTrace, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), GloballyPositioned< T >::position(), Propagator::propagate(), propagator_, sign(), and PV3DBase< T, PVType, FrameType >::z().

◆ propagate() [2/3]

void CachedTrajectory::propagate ( SteppingHelixStateInfo state,
const Cylinder cylinder 
)
private

Definition at line 96 of file CachedTrajectory.cc.

96  {
97  if (const SteppingHelixPropagator* shp = dynamic_cast<const SteppingHelixPropagator*>(propagator_)) {
98  try {
99  shp->propagate(state, cylinder, state);
100  } catch (cms::Exception& ex) {
101  edm::LogWarning("TrackAssociator") << "Caught exception " << ex.category() << ": " << ex.explainSelf();
102  edm::LogWarning("TrackAssociator") << "An exception is caught during the track propagation\n"
103  << state.momentum().x() << ", " << state.momentum().y() << ", "
104  << state.momentum().z();
106  }
107  } else {
109  state.getFreeState(fts);
110  TrajectoryStateOnSurface stateOnSurface = propagator_->propagate(fts, cylinder);
111  state = SteppingHelixStateInfo(*(stateOnSurface.freeState()));
112  }
113 }

References cms::Exception::category(), cms::Exception::explainSelf(), TrajectoryStateOnSurface::freeState(), Propagator::propagate(), and propagator_.

◆ propagate() [3/3]

void CachedTrajectory::propagate ( SteppingHelixStateInfo state,
const Plane plane 
)
private

Definition at line 77 of file CachedTrajectory.cc.

77  {
78  if (const SteppingHelixPropagator* shp = dynamic_cast<const SteppingHelixPropagator*>(propagator_)) {
79  try {
80  shp->propagate(state, plane, state);
81  } catch (cms::Exception& ex) {
82  edm::LogWarning("TrackAssociator") << "Caught exception " << ex.category() << ": " << ex.explainSelf();
83  edm::LogWarning("TrackAssociator") << "An exception is caught during the track propagation\n"
84  << state.momentum().x() << ", " << state.momentum().y() << ", "
85  << state.momentum().z();
87  }
88  } else {
90  state.getFreeState(fts);
91  TrajectoryStateOnSurface stateOnSurface = propagator_->propagate(fts, plane);
92  state = SteppingHelixStateInfo(*(stateOnSurface.freeState()));
93  }
94 }

References cms::Exception::category(), cms::Exception::explainSelf(), TrajectoryStateOnSurface::freeState(), Propagator::propagate(), and propagator_.

Referenced by getTrajectory(), and propagateForward().

◆ propagateAll()

bool CachedTrajectory::propagateAll ( const SteppingHelixStateInfo initialState)
private

propagate through the whole detector, returns true if successful

Definition at line 115 of file CachedTrajectory.cc.

115  {
116  if (fullTrajectoryFilled_) {
117  edm::LogWarning("TrackAssociator")
118  << "Reseting all trajectories. Please call reset_trajectory() explicitely to avoid this message";
120  }
121 
122  // TimerStack timers(TimerStack::Disableable);
123 
125  if (propagator_ == nullptr)
126  throw cms::Exception("FatalError") << "Track propagator is not defined\n";
127  SteppingHelixStateInfo currentState(initialState);
128  fullTrajectory_.push_back(currentState);
129 
130  while (currentState.position().perp() < maxRho_ && fabs(currentState.position().z()) < maxZ_) {
131  LogTrace("TrackAssociator") << "[propagateAll] Propagate outward from (rho, r, z, phi) ("
132  << currentState.position().perp() << ", " << currentState.position().mag() << ", "
133  << currentState.position().z() << ", " << currentState.position().phi() << ")";
134  propagateForward(currentState, step_);
135  if (!currentState.isValid()) {
136  LogTrace("TrackAssociator") << "Failed to propagate the track; moving on\n";
137  break;
138  }
139  LogTrace("TrackAssociator") << "\treached (rho, r, z, phi) (" << currentState.position().perp() << ", "
140  << currentState.position().mag() << ", " << currentState.position().z() << ", "
141  << currentState.position().phi() << ")";
142  fullTrajectory_.push_back(currentState);
143  }
144 
145  SteppingHelixStateInfo currentState2(initialState);
146  SteppingHelixStateInfo previousState;
147  while (currentState2.position().perp() > minRho_ || fabs(currentState2.position().z()) > minZ_) {
148  previousState = currentState2;
149  propagateForward(currentState2, -step_);
150  if (!currentState2.isValid()) {
151  LogTrace("TrackAssociator") << "Failed to propagate the track; moving on\n";
152  break;
153  }
154  if (previousState.position().perp() - currentState2.position().perp() < 0) {
155  LogTrace("TrackAssociator")
156  << "Error: TrackAssociator has propogated the particle past the point of closest approach to IP" << std::endl;
157  break;
158  }
159  LogTrace("TrackAssociator") << "[propagateAll] Propagated inward from (rho, r, z, phi) ("
160  << previousState.position().perp() << ", " << previousState.position().mag() << ", "
161  << previousState.position().z() << "," << previousState.position().phi() << ") to ("
162  << currentState2.position().perp() << ", " << currentState2.position().mag() << ", "
163  << currentState2.position().z() << ", " << currentState2.position().phi() << ")";
164  fullTrajectory_.push_front(currentState2);
165  }
166 
167  // LogTrace("TrackAssociator") << "fullTrajectory_ has " << fullTrajectory_.size() << " states with (R, z):\n";
168  // for(unsigned int i=0; i<fullTrajectory_.size(); i++) {
169  // LogTrace("TrackAssociator") << "state " << i << ": (" << fullTrajectory_[i].position().perp() << ", "
170  // << fullTrajectory_[i].position().z() << ")\n";
171  // }
172 
173  LogTrace("TrackAssociator") << "Done with the track propagation in the detector. Number of steps: "
174  << fullTrajectory_.size();
175  fullTrajectoryFilled_ = true;
176  return !fullTrajectory_.empty();
177 }

References Exception, fullTrajectory_, fullTrajectoryFilled_, SteppingHelixStateInfo::isValid(), LogTrace, PV3DBase< T, PVType, FrameType >::mag(), maxRho_, maxZ_, minRho_, minZ_, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), SteppingHelixStateInfo::position(), propagateForward(), propagator_, reset_trajectory(), step_, and PV3DBase< T, PVType, FrameType >::z().

Referenced by propagateThoughFromIP().

◆ propagateForward()

void CachedTrajectory::propagateForward ( SteppingHelixStateInfo state,
float  distance 
)
private

Definition at line 58 of file CachedTrajectory.cc.

58  {
59  // defined a normal plane wrt the particle trajectory direction
60  // let's hope that I computed the rotation matrix correctly.
61  GlobalVector vector(state.momentum().unit());
62  float r21 = 0;
63  float r22 = vector.z() / sqrt(1 - pow(vector.x(), 2));
64  float r23 = -vector.y() / sqrt(1 - pow(vector.x(), 2));
65  float r31 = vector.x();
66  float r32 = vector.y();
67  float r33 = vector.z();
68  float r11 = r22 * r33 - r23 * r32;
69  float r12 = r23 * r31;
70  float r13 = -r22 * r31;
71 
72  Surface::RotationType rotation(r11, r12, r13, r21, r22, r23, r31, r32, r33);
75 }

References Plane::build(), distance(), funct::pow(), propagate(), dttmaxenums::r32, idealTransformation::rotation, mathSSE::sqrt(), filterCSVwithJSON::target, and trackerHitRTTI::vector.

Referenced by getTrajectory(), and propagateAll().

◆ reset_trajectory()

void CachedTrajectory::reset_trajectory ( )
private

◆ setMaxDetectorLength()

void CachedTrajectory::setMaxDetectorLength ( float  l = 2200.)
inlineprivate

Definition at line 114 of file CachedTrajectory.h.

114 { maxZ_ = l / 2.; }

References cmsLHEtoEOSManager::l, and maxZ_.

Referenced by CachedTrajectory(), and propagateThoughFromIP().

◆ setMaxDetectorRadius()

void CachedTrajectory::setMaxDetectorRadius ( float  r = 800.)
inlineprivate

Definition at line 113 of file CachedTrajectory.h.

113 { maxRho_ = r; }

References maxRho_, and alignCSCRings::r.

Referenced by CachedTrajectory(), and propagateThoughFromIP().

◆ setMaxHOLength()

void CachedTrajectory::setMaxHOLength ( float  l = 2200.)
inlineprivate

Definition at line 116 of file CachedTrajectory.h.

116 { HOmaxZ_ = l / 2.; }

References HOmaxZ_, and cmsLHEtoEOSManager::l.

◆ setMaxHORadius()

void CachedTrajectory::setMaxHORadius ( float  r = 800.)
inlineprivate

Definition at line 115 of file CachedTrajectory.h.

115 { HOmaxRho_ = r; }

References HOmaxRho_, and alignCSCRings::r.

◆ setMinDetectorLength()

void CachedTrajectory::setMinDetectorLength ( float  l = 0.)
inlineprivate

Definition at line 118 of file CachedTrajectory.h.

118 { minZ_ = l / 2.; }

References cmsLHEtoEOSManager::l, and minZ_.

Referenced by CachedTrajectory(), and propagateThoughFromIP().

◆ setMinDetectorRadius()

void CachedTrajectory::setMinDetectorRadius ( float  r = 0.)
inlineprivate

Definition at line 117 of file CachedTrajectory.h.

117 { minRho_ = r; }

References minRho_, and alignCSCRings::r.

Referenced by CachedTrajectory(), and propagateThoughFromIP().

◆ setPropagationStep()

void CachedTrajectory::setPropagationStep ( float  s = 20.)
inlineprivate

Definition at line 120 of file CachedTrajectory.h.

120 { step_ = s; }

References alignCSCRings::s, and step_.

Referenced by CachedTrajectory(), and propagateThoughFromIP().

◆ setPropagator()

void CachedTrajectory::setPropagator ( const Propagator ptr)
inlineprivate

Definition at line 83 of file CachedTrajectory.h.

83 { propagator_ = ptr; }

References propagator_.

Referenced by propagateThoughFromIP().

◆ setStateAtIP()

void CachedTrajectory::setStateAtIP ( const SteppingHelixStateInfo state)
inlineprivate

Definition at line 84 of file CachedTrajectory.h.

84 { stateAtIP_ = state; }

References stateAtIP_.

Referenced by propagateThoughFromIP().

◆ sign()

static int CachedTrajectory::sign ( float  number)
inlinestaticprotected

Definition at line 124 of file CachedTrajectory.h.

124  {
125  if (number == 0)
126  return 0;
127  if (number > 0)
128  return 1;
129  else
130  return -1;
131  }

References contentValuesFiles::number.

Referenced by propagate().

◆ trajectoryDelta()

std::pair< float, float > CachedTrajectory::trajectoryDelta ( TrajectorType  trajectoryType)
private

calculate trajectory change (Theta,Phi) delta = final - original

Definition at line 256 of file CachedTrajectory.cc.

256  {
257  // MEaning of trajectory change depends on its usage. In most cases we measure
258  // change in a trajectory as difference between final track position and initial
259  // direction. In some cases such as change of trajectory in the muon detector we
260  // might want to compare theta-phi of two points or even find local maximum and
261  // mimimum. In general it's not essential what defenition of the trajectory change
262  // is used since we use these numbers only as a rough estimate on how much wider
263  // we should make the preselection region.
264  std::pair<float, float> result(0, 0);
265  if (!stateAtIP_.isValid()) {
266  edm::LogWarning("TrackAssociator") << "State at IP is not known set. Cannot estimate trajectory change. "
267  << "Trajectory change is not taken into account in matching";
268  return result;
269  }
270  switch (trajectoryType) {
271  case IpToEcal:
272  if (ecalTrajectory_.empty())
273  edm::LogWarning("TrackAssociator") << "ECAL trajector is empty. Cannot estimate trajectory change. "
274  << "Trajectory change is not taken into account in matching";
275  else
276  return delta(stateAtIP_.momentum().theta(),
277  ecalTrajectory_.front().position().theta(),
278  stateAtIP_.momentum().phi(),
279  ecalTrajectory_.front().position().phi());
280  break;
281  case IpToHcal:
282  if (hcalTrajectory_.empty())
283  edm::LogWarning("TrackAssociator") << "HCAL trajector is empty. Cannot estimate trajectory change. "
284  << "Trajectory change is not taken into account in matching";
285  else
286  return delta(stateAtIP_.momentum().theta(),
287  hcalTrajectory_.front().position().theta(),
288  stateAtIP_.momentum().phi(),
289  hcalTrajectory_.front().position().phi());
290  break;
291  case IpToHO:
292  if (hoTrajectory_.empty())
293  edm::LogWarning("TrackAssociator") << "HO trajector is empty. Cannot estimate trajectory change. "
294  << "Trajectory change is not taken into account in matching";
295  else
296  return delta(stateAtIP_.momentum().theta(),
297  hoTrajectory_.front().position().theta(),
298  stateAtIP_.momentum().phi(),
299  hoTrajectory_.front().position().phi());
300  break;
301  case FullTrajectory:
302  if (fullTrajectory_.empty())
303  edm::LogWarning("TrackAssociator") << "Full trajector is empty. Cannot estimate trajectory change. "
304  << "Trajectory change is not taken into account in matching";
305  else
306  return delta(stateAtIP_.momentum().theta(),
307  fullTrajectory_.back().position().theta(),
308  stateAtIP_.momentum().phi(),
309  fullTrajectory_.back().position().phi());
310  break;
311  default:
312  edm::LogWarning("TrackAssociator")
313  << "Unkown or not supported trajector type. Cannot estimate trajectory change. "
314  << "Trajectory change is not taken into account in matching";
315  }
316  return result;
317 }

References delta(), ecalTrajectory_, FullTrajectory, fullTrajectory_, hcalTrajectory_, hoTrajectory_, IpToEcal, IpToHcal, IpToHO, SteppingHelixStateInfo::isValid(), SteppingHelixStateInfo::momentum(), PV3DBase< T, PVType, FrameType >::phi(), mps_fire::result, stateAtIP_, and PV3DBase< T, PVType, FrameType >::theta().

Friends And Related Function Documentation

◆ propagateThoughFromIP

std::vector<SteppingHelixStateInfo> propagateThoughFromIP ( const SteppingHelixStateInfo state,
const Propagator ptr,
const FiducialVolume volume,
int  nsteps,
float  step,
float  minR,
float  minZ,
float  maxR,
float  maxZ 
)
friend

Definition at line 18 of file CachedTrajectory.cc.

26  {
27  CachedTrajectory neckLace;
28  neckLace.setStateAtIP(state);
29  neckLace.reset_trajectory();
30  neckLace.setPropagator(prop);
31  neckLace.setPropagationStep(0.1);
32  neckLace.setMinDetectorRadius(minR);
33  neckLace.setMinDetectorLength(minZ * 2.);
34  neckLace.setMaxDetectorRadius(maxR);
35  neckLace.setMaxDetectorLength(maxZ * 2.);
36 
37  // Propagate track
38  bool isPropagationSuccessful = neckLace.propagateAll(state);
39 
40  if (!isPropagationSuccessful)
41  return std::vector<SteppingHelixStateInfo>();
42 
43  std::vector<SteppingHelixStateInfo> complicatePoints;
44  neckLace.getTrajectory(complicatePoints, volume, nsteps);
45 
46  return complicatePoints;
47 }

◆ TrackDetectorAssociator

friend class TrackDetectorAssociator
friend

Definition at line 52 of file CachedTrajectory.h.

Member Data Documentation

◆ ecalTrajectory_

std::vector<SteppingHelixStateInfo> CachedTrajectory::ecalTrajectory_
protected

◆ fullTrajectory_

std::deque<SteppingHelixStateInfo> CachedTrajectory::fullTrajectory_
protected

◆ fullTrajectoryFilled_

bool CachedTrajectory::fullTrajectoryFilled_
protected

Definition at line 152 of file CachedTrajectory.h.

Referenced by getTrajectory(), propagateAll(), and reset_trajectory().

◆ hcalTrajectory_

std::vector<SteppingHelixStateInfo> CachedTrajectory::hcalTrajectory_
protected

◆ HOmaxRho_

float CachedTrajectory::HOmaxRho_
protected

Definition at line 158 of file CachedTrajectory.h.

Referenced by setMaxHORadius().

◆ HOmaxZ_

float CachedTrajectory::HOmaxZ_
protected

Definition at line 159 of file CachedTrajectory.h.

Referenced by setMaxHOLength().

◆ hoTrajectory_

std::vector<SteppingHelixStateInfo> CachedTrajectory::hoTrajectory_
protected

◆ maxRho_

float CachedTrajectory::maxRho_
protected

Definition at line 156 of file CachedTrajectory.h.

Referenced by propagateAll(), and setMaxDetectorRadius().

◆ maxZ_

float CachedTrajectory::maxZ_
protected

Definition at line 157 of file CachedTrajectory.h.

Referenced by propagateAll(), and setMaxDetectorLength().

◆ minRho_

float CachedTrajectory::minRho_
protected

Definition at line 160 of file CachedTrajectory.h.

Referenced by propagateAll(), and setMinDetectorRadius().

◆ minZ_

float CachedTrajectory::minZ_
protected

Definition at line 161 of file CachedTrajectory.h.

Referenced by propagateAll(), and setMinDetectorLength().

◆ preshowerTrajectory_

std::vector<SteppingHelixStateInfo> CachedTrajectory::preshowerTrajectory_
protected

◆ propagator_

const Propagator* CachedTrajectory::propagator_
protected

Definition at line 154 of file CachedTrajectory.h.

Referenced by propagate(), propagateAll(), and setPropagator().

◆ stateAtIP_

SteppingHelixStateInfo CachedTrajectory::stateAtIP_
protected

Definition at line 150 of file CachedTrajectory.h.

Referenced by setStateAtIP(), and trajectoryDelta().

◆ step_

float CachedTrajectory::step_
protected

Definition at line 162 of file CachedTrajectory.h.

Referenced by getPropagationStep(), propagateAll(), and setPropagationStep().

◆ wideEcalTrajectory_

std::vector<GlobalPoint> CachedTrajectory::wideEcalTrajectory_
protected

Definition at line 147 of file CachedTrajectory.h.

Referenced by getWideTrajectory(), and reset_trajectory().

◆ wideHcalTrajectory_

std::vector<GlobalPoint> CachedTrajectory::wideHcalTrajectory_
protected

Definition at line 148 of file CachedTrajectory.h.

Referenced by getWideTrajectory(), and reset_trajectory().

◆ wideHOTrajectory_

std::vector<GlobalPoint> CachedTrajectory::wideHOTrajectory_
protected

Definition at line 149 of file CachedTrajectory.h.

Referenced by getWideTrajectory(), and reset_trajectory().

Vector3DBase
Definition: Vector3DBase.h:8
Point2DBase< float, LocalTag >
SteppingHelixStateInfo::momentum
GlobalVector momentum() const
Definition: SteppingHelixStateInfo.h:60
TkRotation< float >
mps_fire.i
i
Definition: mps_fire.py:428
Reference_intrackfit_cff.barrel
list barrel
Definition: Reference_intrackfit_cff.py:37
SteppingHelixStateInfo::position
GlobalPoint position() const
Definition: SteppingHelixStateInfo.h:59
HLT_FULL_cff.minR
minR
Definition: HLT_FULL_cff.py:95122
CachedTrajectory::HOmaxZ_
float HOmaxZ_
Definition: CachedTrajectory.h:159
makePileupJSON.denom
denom
Definition: makePileupJSON.py:147
step
step
Definition: StallMonitor.cc:94
LocalError::xy
float xy() const
Definition: LocalError.h:23
CachedTrajectory::setMaxDetectorRadius
void setMaxDetectorRadius(float r=800.)
Definition: CachedTrajectory.h:113
CachedTrajectory::setMinDetectorLength
void setMinDetectorLength(float l=0.)
Definition: CachedTrajectory.h:118
CachedTrajectory::wideHOTrajectory_
std::vector< GlobalPoint > wideHOTrajectory_
Definition: CachedTrajectory.h:149
CachedTrajectory::HO
Definition: CachedTrajectory.h:65
CachedTrajectory::setMinDetectorRadius
void setMinDetectorRadius(float r=0.)
Definition: CachedTrajectory.h:117
PV3DBase::theta
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
FiducialVolume::minR
double minR(bool withTolerance=true) const
Definition: FiducialVolume.h:37
CachedTrajectory::FullTrajectory
Definition: CachedTrajectory.h:64
CachedTrajectory::maxRho_
float maxRho_
Definition: CachedTrajectory.h:156
CachedTrajectory::getTrajectory
void getTrajectory(std::vector< SteppingHelixStateInfo > &, const FiducialVolume &, int steps=4)
Definition: CachedTrajectory.cc:334
ReferenceCountingPointer< Cylinder >
makeMuonMisalignmentScenario.endcap
endcap
Definition: makeMuonMisalignmentScenario.py:320
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
GloballyPositioned< float >::RotationType
TkRotation< float > RotationType
Definition: GloballyPositioned.h:22
CachedTrajectory::IpToEcal
Definition: CachedTrajectory.h:64
CachedTrajectory::setPropagationStep
void setPropagationStep(float s=20.)
Definition: CachedTrajectory.h:120
CachedTrajectory::setMaxDetectorLength
void setMaxDetectorLength(float l=2200.)
Definition: CachedTrajectory.h:114
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
CosmicsPD_Skims.maxZ
maxZ
Definition: CosmicsPD_Skims.py:136
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
CachedTrajectory::minRho_
float minRho_
Definition: CachedTrajectory.h:160
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
alignCSCRings.s
s
Definition: alignCSCRings.py:92
SteppingHelixPropagator
Definition: SteppingHelixPropagator.h:36
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
TrajectoryStateOnSurface::freeState
FreeTrajectoryState const * freeState(bool withErrors=true) const
Definition: TrajectoryStateOnSurface.h:58
contentValuesFiles.number
number
Definition: contentValuesFiles.py:53
CachedTrajectory::propagateAll
bool propagateAll(const SteppingHelixStateInfo &initialState)
propagate through the whole detector, returns true if successful
Definition: CachedTrajectory.cc:115
LocalError::xx
float xx() const
Definition: LocalError.h:22
CachedTrajectory::setPropagator
void setPropagator(const Propagator *ptr)
Definition: CachedTrajectory.h:83
FiducialVolume::maxR
double maxR(bool withTolerance=true) const
Definition: FiducialVolume.h:43
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
CachedTrajectory::propagate
void propagate(SteppingHelixStateInfo &state, const Plane &plane)
Definition: CachedTrajectory.cc:77
DDAxes::z
Point3DBase< float, GlobalTag >
CachedTrajectory::propagator_
const Propagator * propagator_
Definition: CachedTrajectory.h:154
geometryCSVtoXML.xy
xy
Definition: geometryCSVtoXML.py:19
LocalError::scale
LocalError scale(float s) const
Definition: LocalError.h:31
SteppingHelixStateInfo
Definition: SteppingHelixStateInfo.h:27
Plane::build
static PlanePointer build(Args &&... args)
Definition: Plane.h:33
Geom::Phi::phi
T1 phi() const
Definition: Phi.h:78
idealTransformation.rotation
dictionary rotation
Definition: idealTransformation.py:1
CachedTrajectory::wideHcalTrajectory_
std::vector< GlobalPoint > wideHcalTrajectory_
Definition: CachedTrajectory.h:148
geometryCSVtoXML.yy
yy
Definition: geometryCSVtoXML.py:19
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
FiducialVolume::minZ
double minZ(bool withTolerance=true) const
Definition: FiducialVolume.h:49
LocalError
Definition: LocalError.h:12
CachedTrajectory::fullTrajectory_
std::deque< SteppingHelixStateInfo > fullTrajectory_
Definition: CachedTrajectory.h:142
CachedTrajectory::sign
static int sign(float number)
Definition: CachedTrajectory.h:124
LocalError::rotate
LocalError rotate(float x, float y) const
Return a new LocalError, rotated by an angle defined by the direction (x,y)
Definition: LocalError.h:37
CachedTrajectory::IpToHO
Definition: CachedTrajectory.h:64
CachedTrajectory::maxZ_
float maxZ_
Definition: CachedTrajectory.h:157
Plane::localZ
float localZ(const GlobalPoint &gp) const
Definition: Plane.h:45
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
PV3DBase::eta
T eta() const
Definition: PV3DBase.h:73
CachedTrajectory::wideEcalTrajectory_
std::vector< GlobalPoint > wideEcalTrajectory_
Definition: CachedTrajectory.h:147
Cylinder::build
static CylinderPointer build(const PositionType &pos, const RotationType &rot, Scalar radius, Bounds *bounds=nullptr)
Definition: Cylinder.h:45
createfilelist.int
int
Definition: createfilelist.py:10
FiducialVolume::isValid
bool isValid() const
check whether the volume is properly defined
Definition: FiducialVolume.cc:19
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
CachedTrajectory::hoTrajectory_
std::vector< SteppingHelixStateInfo > hoTrajectory_
Definition: CachedTrajectory.h:145
Propagator::propagate
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:50
CachedTrajectory
Definition: CachedTrajectory.h:44
SteppingHelixStateInfo::isValid
bool isValid() const
Definition: SteppingHelixStateInfo.h:65
CachedTrajectory::minZ_
float minZ_
Definition: CachedTrajectory.h:161
CachedTrajectory::delta
std::pair< float, float > delta(const double &theta1, const double &theta2, const double &phi1, const double &phi2)
Definition: CachedTrajectory.cc:319
cmsLHEtoEOSManager.l
l
Definition: cmsLHEtoEOSManager.py:204
CachedTrajectory::Ecal
Definition: CachedTrajectory.h:65
alignCSCRings.r
r
Definition: alignCSCRings.py:93
FiducialVolume::maxZ
double maxZ(bool withTolerance=true) const
Definition: FiducialVolume.h:55
DDAxes::phi
PV3DBase::mag
T mag() const
Definition: PV3DBase.h:64
GloballyPositioned::position
const PositionType & position() const
Definition: GloballyPositioned.h:36
CachedTrajectory::propagateForward
void propagateForward(SteppingHelixStateInfo &state, float distance)
Definition: CachedTrajectory.cc:58
CachedTrajectory::IpToHcal
Definition: CachedTrajectory.h:64
RunInfoPI::state
state
Definition: RunInfoPayloadInspectoHelper.h:16
GloballyPositioned< float >::PositionType
Point3DBase< float, GlobalTag > PositionType
Definition: GloballyPositioned.h:21
FreeTrajectoryState
Definition: FreeTrajectoryState.h:27
dttmaxenums::r32
Definition: DTTMax.h:28
CachedTrajectory::fullTrajectoryFilled_
bool fullTrajectoryFilled_
Definition: CachedTrajectory.h:152
CachedTrajectory::setStateAtIP
void setStateAtIP(const SteppingHelixStateInfo &state)
Definition: CachedTrajectory.h:84
CachedTrajectory::step_
float step_
Definition: CachedTrajectory.h:162
Exception
Definition: hltDiff.cc:245
Plane
Definition: Plane.h:16
CachedTrajectory::preshowerTrajectory_
std::vector< SteppingHelixStateInfo > preshowerTrajectory_
Definition: CachedTrajectory.h:146
cms::Exception::explainSelf
virtual std::string explainSelf() const
Definition: Exception.cc:108
CachedTrajectory::stateAtIP_
SteppingHelixStateInfo stateAtIP_
Definition: CachedTrajectory.h:150
CachedTrajectory::distance
float distance(const Plane *plane, int index)
Definition: CachedTrajectory.h:136
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
filterCSVwithJSON.target
target
Definition: filterCSVwithJSON.py:32
customisers.steps
steps
Definition: customisers.py:40
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
CachedTrajectory::ecalTrajectory_
std::vector< SteppingHelixStateInfo > ecalTrajectory_
Definition: CachedTrajectory.h:143
M_PI_2
#define M_PI_2
Definition: CSCGattiFunction.cc:6
mps_fire.result
result
Definition: mps_fire.py:311
cms::Exception
Definition: Exception.h:70
LogTrace
#define LogTrace(id)
Definition: MessageLogger.h:234
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition: HGC3DClusterGenMatchSelector_cfi.py:7
Local2DPoint
Point2DBase< float, LocalTag > Local2DPoint
Definition: LocalPoint.h:8
cms::Exception::category
std::string const & category() const
Definition: Exception.cc:143
edm::Log
Definition: MessageLogger.h:70
CachedTrajectory::reset_trajectory
void reset_trajectory()
Definition: CachedTrajectory.cc:431
CachedTrajectory::hcalTrajectory_
std::vector< SteppingHelixStateInfo > hcalTrajectory_
Definition: CachedTrajectory.h:144
LocalError::yy
float yy() const
Definition: LocalError.h:24
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
geometryCSVtoXML.xx
xx
Definition: geometryCSVtoXML.py:19
CachedTrajectory::HOmaxRho_
float HOmaxRho_
Definition: CachedTrajectory.h:158
CosmicGenFilterHelix_cff.minZ
minZ
i.e.
Definition: CosmicGenFilterHelix_cff.py:9
CachedTrajectory::Hcal
Definition: CachedTrajectory.h:65