CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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, edm::ConsumesCollector iC)
 

Private Attributes

TkClonerImpl theHitCloner
 
const int theNumberMeasurementsForFit
 
const PropagatorthePropagatorAlong
 
const edm::ESGetToken
< Propagator,
TrackingComponentsRecord
thePropagatorAlongToken
 
const PropagatorthePropagatorOpposite
 
const edm::ESGetToken
< Propagator,
TrackingComponentsRecord
thePropagatorOppositeToken
 

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 24 of file TransientInitialStateEstimator.h.

Member Typedef Documentation

Definition at line 26 of file TransientInitialStateEstimator.h.

Constructor & Destructor Documentation

TransientInitialStateEstimator::TransientInitialStateEstimator ( const edm::ParameterSet conf,
edm::ConsumesCollector  iC 
)

Definition at line 20 of file TransientInitialStateEstimator.cc.

22  iC.esConsumes(edm::ESInputTag("", conf.getParameter<std::string>("propagatorAlongTISE")))),
24  iC.esConsumes(edm::ESInputTag("", conf.getParameter<std::string>("propagatorOppositeTISE")))),
25  thePropagatorAlong(nullptr),
26  thePropagatorOpposite(nullptr),
27  theNumberMeasurementsForFit(conf.getParameter<int>("numberMeasurementsForFit")) {}
const edm::ESGetToken< Propagator, TrackingComponentsRecord > thePropagatorAlongToken
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::ESGetToken< Propagator, TrackingComponentsRecord > thePropagatorOppositeToken

Member Function Documentation

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

Definition at line 35 of file TransientInitialStateEstimator.cc.

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

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

Member Data Documentation

TkClonerImpl TransientInitialStateEstimator::theHitCloner
private

Definition at line 38 of file TransientInitialStateEstimator.h.

Referenced by innerState(), and setEventSetup().

const int TransientInitialStateEstimator::theNumberMeasurementsForFit
private

Definition at line 39 of file TransientInitialStateEstimator.h.

Referenced by innerState().

const Propagator* TransientInitialStateEstimator::thePropagatorAlong
private

Definition at line 36 of file TransientInitialStateEstimator.h.

Referenced by innerState(), and setEventSetup().

const edm::ESGetToken<Propagator, TrackingComponentsRecord> TransientInitialStateEstimator::thePropagatorAlongToken
private

Definition at line 34 of file TransientInitialStateEstimator.h.

Referenced by setEventSetup().

const Propagator* TransientInitialStateEstimator::thePropagatorOpposite
private

Definition at line 37 of file TransientInitialStateEstimator.h.

Referenced by setEventSetup().

const edm::ESGetToken<Propagator, TrackingComponentsRecord> TransientInitialStateEstimator::thePropagatorOppositeToken
private

Definition at line 35 of file TransientInitialStateEstimator.h.

Referenced by setEventSetup().