CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
TransientInitialStateEstimator Class Reference

#include <TransientInitialStateEstimator.h>

Public Types

typedef TrajectoryStateOnSurface TSOS
 

Public Member Functions

std::pair< TrajectoryStateOnSurface, const GeomDet * > innerState (const Trajectory &traj, bool doBackFit=true) const
 
void setEventSetup (const edm::EventSetup &es, const TkClonerImpl &hc)
 
 TransientInitialStateEstimator (const edm::ParameterSet &conf)
 

Private Attributes

TkClonerImpl theHitCloner
 
const int theNumberMeasurementsForFit
 
const PropagatorthePropagatorAlong
 
const std::string thePropagatorAlongName
 
const PropagatorthePropagatorOpposite
 
const std::string thePropagatorOppositeName
 

Detailed Description

Computes the trajectory state to be used as a starting state for the track fit from the vector of hits. The parameters of this state are close to the final fit parameters. The error matrix is enlarged in order not to bias the track fit.

Definition at line 23 of file TransientInitialStateEstimator.h.

Member Typedef Documentation

Definition at line 25 of file TransientInitialStateEstimator.h.

Constructor & Destructor Documentation

TransientInitialStateEstimator::TransientInitialStateEstimator ( const edm::ParameterSet conf)

Definition at line 21 of file TransientInitialStateEstimator.cc.

22  : thePropagatorAlongName(conf.getParameter<std::string>("propagatorAlongTISE")),
23  thePropagatorOppositeName(conf.getParameter<std::string>("propagatorOppositeTISE")),
24  thePropagatorAlong(nullptr),
25  thePropagatorOpposite(nullptr),
26  theNumberMeasurementsForFit(conf.getParameter<int>("numberMeasurementsForFit")) {}
T getParameter(std::string const &) const

Member Function Documentation

std::pair< TrajectoryStateOnSurface, const GeomDet * > TransientInitialStateEstimator::innerState ( const Trajectory traj,
bool  doBackFit = true 
) const

Definition at line 40 of file TransientInitialStateEstimator.cc.

References alongMomentum, Chi2MeasurementEstimator_cfi::Chi2MeasurementEstimator, Trajectory::direction(), Trajectory::firstMeasurement(), TrajectoryMeasurement::forwardPredictedState(), mps_fire::i, TrajectoryStateOnSurface::isValid(), Trajectory::isValid(), Trajectory::lastMeasurement(), LogDebug, TrajectoryFitter::looper, Propagator::magneticField(), Trajectory::measurements(), eostools::move(), Trajectory::nLoops(), oppositeToMomentum, rpcPointValidation_cfi::recHit, TrajectoryMeasurement::recHit(), TrajectoryStateOnSurface::rescaleError(), TrajectoryFitter::standard, theHitCloner, theNumberMeasurementsForFit, and thePropagatorAlong.

41  {
42  if (!doBackFit && traj.firstMeasurement().forwardPredictedState().isValid()) {
43  LogDebug("TransientInitialStateEstimator")
44  << "a backward fit will not be done. assuming that the state on first measurement is OK";
45  TSOS firstStateFromForward = traj.firstMeasurement().forwardPredictedState();
46  firstStateFromForward.rescaleError(100.);
47  return std::pair<TrajectoryStateOnSurface, const GeomDet*>(std::move(firstStateFromForward),
48  traj.firstMeasurement().recHit()->det());
49  }
50  if (!doBackFit) {
51  LogDebug("TransientInitialStateEstimator")
52  << "told to not do a back fit, but the forward state of the first measurement is not valid. doing a back fit.";
53  }
54 
55  int nMeas = traj.measurements().size();
56  int lastFitted = theNumberMeasurementsForFit >= 0 ? theNumberMeasurementsForFit : nMeas;
57  if (nMeas - 1 < lastFitted)
58  lastFitted = nMeas - 1;
59 
60  auto const& measvec = traj.measurements();
62 
63  bool foundLast = false;
64  int actualLast = -99;
65 
66  for (int i = lastFitted; i >= 0; i--) {
67  if (measvec[i].recHit()->det()) {
68  if (!foundLast) {
69  actualLast = i;
70  foundLast = true;
71  }
72  firstHits.push_back(measvec[i].recHit());
73  }
74  }
75  TSOS startingState = measvec[actualLast].updatedState();
76  startingState.rescaleError(100.);
77 
78  // avoid cloning...
79  KFUpdator const aKFUpdator;
80  Chi2MeasurementEstimator const aChi2MeasurementEstimator(100., 3);
81  KFTrajectoryFitter backFitter(
82  thePropagatorAlong, &aKFUpdator, &aChi2MeasurementEstimator, firstHits.size(), nullptr, &theHitCloner);
83 
85 
86  // only direction matters in this contest
88 
89  Trajectory&& fitres = backFitter.fitOne(
90  fakeSeed, firstHits, startingState, traj.nLoops() > 0 ? TrajectoryFitter::looper : TrajectoryFitter::standard);
91 
92  LogDebug("TransientInitialStateEstimator")
93  << "using a backward fit of :" << firstHits.size() << " hits, starting from:\n"
94  << startingState << " to get the estimate of the initial state of the track.";
95 
96  if (!fitres.isValid()) {
97  LogDebug("TransientInitialStateEstimator") << "FitTester: first hits fit failed!";
98  return std::pair<TrajectoryStateOnSurface, const GeomDet*>();
99  }
100 
101  TrajectoryMeasurement const& firstMeas = fitres.lastMeasurement();
102 
103  // magnetic field can be different!
104  TSOS firstState(firstMeas.updatedState().localParameters(),
105  firstMeas.updatedState().localError(),
106  firstMeas.updatedState().surface(),
108 
109  // TSOS & firstState = const_cast<TSOS&>(firstMeas.updatedState());
110 
111  // this fails in case of zero field? (for sure for beamhalo reconstruction)
112  // assert(thePropagatorAlong->magneticField()==firstState.magneticField());
113 
114  firstState.rescaleError(100.);
115 
116  LogDebug("TransientInitialStateEstimator")
117  << "the initial state is found to be:\n:" << firstState
118  << "\n it's field pointer is: " << firstState.magneticField()
119  << "\n the pointer from the state of the back fit was: " << firstMeas.updatedState().magneticField();
120 
121  return std::pair<TrajectoryStateOnSurface, const GeomDet*>(std::move(firstState), firstMeas.recHit()->det());
122 }
#define LogDebug(id)
ConstRecHitPointer const & recHit() const
PropagationDirection
signed char nLoops() const
Definition: Trajectory.h:329
PropagationDirection const & direction() const
Definition: Trajectory.cc:133
DataContainer const & measurements() const
Definition: Trajectory.h:178
virtual const MagneticField * magneticField() const =0
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:150
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
std::vector< ConstRecHitPointer > ConstRecHitContainer
bool isValid() const
Definition: Trajectory.h:257
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:166
def move(src, dest)
Definition: eostools.py:511
void TransientInitialStateEstimator::setEventSetup ( const edm::EventSetup es,
const TkClonerImpl hc 
)

Definition at line 28 of file TransientInitialStateEstimator.cc.

References edm::EventSetup::get(), AnalysisDataFormats_SUSYBSMObjects::hc, edm::ESHandle< T >::product(), theHitCloner, thePropagatorAlong, thePropagatorAlongName, thePropagatorOpposite, and thePropagatorOppositeName.

Member Data Documentation

TkClonerImpl TransientInitialStateEstimator::theHitCloner
private

Definition at line 37 of file TransientInitialStateEstimator.h.

Referenced by innerState(), and setEventSetup().

const int TransientInitialStateEstimator::theNumberMeasurementsForFit
private

Definition at line 38 of file TransientInitialStateEstimator.h.

Referenced by innerState().

const Propagator* TransientInitialStateEstimator::thePropagatorAlong
private

Definition at line 35 of file TransientInitialStateEstimator.h.

Referenced by innerState(), and setEventSetup().

const std::string TransientInitialStateEstimator::thePropagatorAlongName
private

Definition at line 33 of file TransientInitialStateEstimator.h.

Referenced by setEventSetup().

const Propagator* TransientInitialStateEstimator::thePropagatorOpposite
private

Definition at line 36 of file TransientInitialStateEstimator.h.

Referenced by setEventSetup().

const std::string TransientInitialStateEstimator::thePropagatorOppositeName
private

Definition at line 34 of file TransientInitialStateEstimator.h.

Referenced by setEventSetup().