CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

TransientInitialStateEstimator Class Reference

#include <TransientInitialStateEstimator.h>

List of all members.

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)
 Call this at each event until this object will come from the EventSetup as it should.
 TransientInitialStateEstimator (const edm::EventSetup &es, const edm::ParameterSet &conf)

Private Attributes

int theNumberMeasurementsForFit
edm::ESHandle< PropagatorthePropagatorAlong
std::string thePropagatorAlongName
edm::ESHandle< PropagatorthePropagatorOpposite
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 20 of file TransientInitialStateEstimator.h.


Member Typedef Documentation

Definition at line 23 of file TransientInitialStateEstimator.h.


Constructor & Destructor Documentation

TransientInitialStateEstimator::TransientInitialStateEstimator ( const edm::EventSetup es,
const edm::ParameterSet conf 
)

Definition at line 21 of file TransientInitialStateEstimator.cc.

References edm::EventSetup::get(), edm::ParameterSet::getParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

{
  thePropagatorAlongName    = conf.getParameter<std::string>("propagatorAlongTISE");   
  thePropagatorOppositeName = conf.getParameter<std::string>("propagatorOppositeTISE");   
  theNumberMeasurementsForFit = conf.getParameter<int32_t>("numberMeasurementsForFit");   


  // let's avoid breaking compatibility now
  es.get<TrackingComponentsRecord>().get(thePropagatorAlongName,thePropagatorAlong);
  es.get<TrackingComponentsRecord>().get(thePropagatorOppositeName,thePropagatorOpposite);
}

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, Trajectory::direction(), Trajectory::firstMeasurement(), TrajectoryMeasurement::forwardPredictedState(), i, TrajectoryStateOnSurface::isValid(), TrajectoryStateOnSurface::localError(), TrajectoryStateOnSurface::localParameters(), LogDebug, TrajectoryStateOnSurface::magneticField(), Trajectory::measurements(), oppositeToMomentum, TrajectoryMeasurement::recHit(), TrajectoryStateOnSurface::rescaleError(), TrajectoryStateOnSurface::surface(), and TrajectoryMeasurement::updatedState().

Referenced by cms::CkfTrackCandidateMakerBase::produceBase(), InOutConversionTrackFinder::tracks(), and OutInConversionTrackFinder::tracks().

{
  if (!doBackFit && traj.firstMeasurement().forwardPredictedState().isValid()){
    LogDebug("TransientInitialStateEstimator")
      <<"a backward fit will not be done. assuming that the state on first measurement is OK";
    TSOS firstStateFromForward = traj.firstMeasurement().forwardPredictedState();
    firstStateFromForward.rescaleError(100.);    
    return std::pair<TrajectoryStateOnSurface, const GeomDet*>( firstStateFromForward, 
                                                                traj.firstMeasurement().recHit()->det());
  }
  if (!doBackFit){
    LogDebug("TransientInitialStateEstimator")
      <<"told to not do a back fit, but the forward state of the first measurement is not valid. doing a back fit.";
  }

  int nMeas = traj.measurements().size();
  int lastFitted = theNumberMeasurementsForFit >=0 ? theNumberMeasurementsForFit : nMeas; 
  if (nMeas-1 < lastFitted) lastFitted = nMeas-1;

  std::vector<TrajectoryMeasurement> measvec = traj.measurements();
  TransientTrackingRecHit::ConstRecHitContainer firstHits;

  bool foundLast = false;
  int actualLast = -99;

  for (int i=lastFitted; i >= 0; i--) {
    if(measvec[i].recHit()->det()){
      if(!foundLast){
        actualLast = i; 
        foundLast = true;
      }
      firstHits.push_back( measvec[i].recHit());
    }
  }
  TSOS startingState = measvec[actualLast].updatedState();
  startingState.rescaleError(100.);

  // avoid cloning...
  KFUpdator const aKFUpdator;
  Chi2MeasurementEstimator const aChi2MeasurementEstimator( 100., 3);
  KFTrajectoryFitter backFitter( thePropagatorAlong.product(),
                                 &aKFUpdator,
                                 &aChi2MeasurementEstimator,
                                 firstHits.size());

  PropagationDirection backFitDirection = traj.direction() == alongMomentum ? oppositeToMomentum: alongMomentum;

  // only direction matters in this contest
  TrajectorySeed fakeSeed = TrajectorySeed(PTrajectoryStateOnDet() , 
                                           edm::OwnVector<TrackingRecHit>(),
                                           backFitDirection);

  vector<Trajectory> fitres = backFitter.fit( fakeSeed, firstHits, startingState);
  
  LogDebug("TransientInitialStateEstimator")
    <<"using a backward fit of :"<<firstHits.size()<<" hits, starting from:\n"<<startingState
    <<" to get the estimate of the initial state of the track.";

  if (fitres.size() != 1) {
        LogDebug("TransientInitialStateEstimator")
          << "FitTester: first hits fit failed!";
    return std::pair<TrajectoryStateOnSurface, const GeomDet*>();
  }

  TrajectoryMeasurement firstMeas = fitres[0].lastMeasurement();
  TSOS firstState(firstMeas.updatedState().localParameters(),
                  firstMeas.updatedState().localError(),
                  firstMeas.updatedState().surface(),
                  thePropagatorAlong->magneticField());
  // I couldn't do: 
  //TSOS firstState = firstMeas.updatedState();
  // why????


  firstState.rescaleError(100.);

  LogDebug("TransientInitialStateEstimator")
    <<"the initial state is found to be:\n:"<<firstState
    <<"\n it's field pointer is: "<<firstState.magneticField()
    <<"\n the pointer from the state of the back fit was: "<<firstMeas.updatedState().magneticField();


  return std::pair<TrajectoryStateOnSurface, const GeomDet*>( firstState, 
                                                              firstMeas.recHit()->det());
}
void TransientInitialStateEstimator::setEventSetup ( const edm::EventSetup es)

Member Data Documentation

Definition at line 38 of file TransientInitialStateEstimator.h.

Definition at line 36 of file TransientInitialStateEstimator.h.

Definition at line 34 of file TransientInitialStateEstimator.h.

Definition at line 37 of file TransientInitialStateEstimator.h.

Definition at line 35 of file TransientInitialStateEstimator.h.