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

Member Typedef Documentation

Definition at line 24 of file TransientInitialStateEstimator.h.

Constructor & Destructor Documentation

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

Definition at line 21 of file TransientInitialStateEstimator.cc.

21  :
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"))
27 {}
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 42 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.

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

Definition at line 29 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 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 std::string TransientInitialStateEstimator::thePropagatorAlongName
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 std::string TransientInitialStateEstimator::thePropagatorOppositeName
private

Definition at line 35 of file TransientInitialStateEstimator.h.

Referenced by setEventSetup().