CMS 3D CMS Logo

KalmanAlignmentTrackRefitter Class Reference

This class serves the very specific needs of the KalmanAlignmentAlgorithm. More...

#include <Alignment/KalmanAlignmentAlgorithm/interface/KalmanAlignmentTrackRefitter.h>

Inheritance diagram for KalmanAlignmentTrackRefitter:

TrackProducerBase< reco::Track >

List of all members.

Public Types

typedef std::vector
< KalmanAlignmentSetup * > 
AlignmentSetupCollection
typedef
AlignmentAlgorithmBase::ConstTrajTrackPair 
ConstTrajTrackPair
typedef
AlignmentAlgorithmBase::ConstTrajTrackPairCollection 
ConstTrajTrackPairCollection
typedef edm::OwnVector
< TrackingRecHit
RecHitContainer
typedef
KalmanAlignmentSetup::SortingDirection 
SortingDirection
typedef std::vector< TrackletPtrTrackletCollection
typedef
KalmanAlignmentTracklet::TrackletPtr 
TrackletPtr
typedef
KalmanAlignmentTracklet::TrajTrackPairCollection 
TrajTrackPairCollection

Public Member Functions

 KalmanAlignmentTrackRefitter (const edm::ParameterSet &config, AlignableNavigator *navigator)
 Constructor.
virtual void produce (edm::Event &, const edm::EventSetup &)
 Dummy implementation, due to inheritance from TrackProducerBase.
TrackletCollection refitTracks (const edm::EventSetup &eventSetup, const AlignmentSetupCollection &algoSetups, const ConstTrajTrackPairCollection &tracks)
 ~KalmanAlignmentTrackRefitter (void)
 Destructor.

Private Member Functions

void debugTrackData (const std::string identifier, const Trajectory *traj, const reco::Track *track)
TrajTrackPairCollection refitSingleTracklet (const TrackingGeometry *geometry, const MagneticField *magneticField, const TrajectoryFitter *fitter, const Propagator *propagator, const TransientTrackingRecHitBuilder *recHitBuilder, const reco::TransientTrack &originalTrack, RecHitContainer &recHits, const SortingDirection &sortingDir, bool useExternalEstimate, bool reuseMomentumEstimate)
void sortRecHits (RecHitContainer &hits, const TransientTrackingRecHitBuilder *builder, const SortingDirection &sortingDir) const

Private Attributes

bool theDebugFlag
AlignableNavigatortheNavigator
TrackProducerAlgorithm
< reco::Track
theRefitterAlgo


Detailed Description

This class serves the very specific needs of the KalmanAlignmentAlgorithm.

Tracks are partially refitted to 'tracklets' using the current estimate on the alignment (see class CurrentAlignmentKFUpdator. These tracklets are either used to compute an exteranal estimate for other tracklets or are handed to the alignment algorithm for further processing. If a tracklet is used as an external prediction or for further processing is defined via the configuration file. NOTE: The trajectory measurements of the tracklets are always ordered along the direction of the momentum!

Definition at line 29 of file KalmanAlignmentTrackRefitter.h.


Member Typedef Documentation

typedef std::vector< KalmanAlignmentSetup* > KalmanAlignmentTrackRefitter::AlignmentSetupCollection

Definition at line 34 of file KalmanAlignmentTrackRefitter.h.

typedef AlignmentAlgorithmBase::ConstTrajTrackPair KalmanAlignmentTrackRefitter::ConstTrajTrackPair

Definition at line 43 of file KalmanAlignmentTrackRefitter.h.

typedef AlignmentAlgorithmBase::ConstTrajTrackPairCollection KalmanAlignmentTrackRefitter::ConstTrajTrackPairCollection

Definition at line 44 of file KalmanAlignmentTrackRefitter.h.

typedef edm::OwnVector< TrackingRecHit > KalmanAlignmentTrackRefitter::RecHitContainer

Definition at line 37 of file KalmanAlignmentTrackRefitter.h.

typedef KalmanAlignmentSetup::SortingDirection KalmanAlignmentTrackRefitter::SortingDirection

Definition at line 35 of file KalmanAlignmentTrackRefitter.h.

typedef std::vector< TrackletPtr > KalmanAlignmentTrackRefitter::TrackletCollection

Definition at line 41 of file KalmanAlignmentTrackRefitter.h.

typedef KalmanAlignmentTracklet::TrackletPtr KalmanAlignmentTrackRefitter::TrackletPtr

Definition at line 40 of file KalmanAlignmentTrackRefitter.h.

typedef KalmanAlignmentTracklet::TrajTrackPairCollection KalmanAlignmentTrackRefitter::TrajTrackPairCollection

Definition at line 39 of file KalmanAlignmentTrackRefitter.h.


Constructor & Destructor Documentation

KalmanAlignmentTrackRefitter::KalmanAlignmentTrackRefitter ( const edm::ParameterSet config,
AlignableNavigator navigator 
)

Constructor.

Definition at line 23 of file KalmanAlignmentTrackRefitter.cc.

References edm::ParameterSet::getParameter(), TrackProducerBase< T >::setConf(), and TrackProducerBase< T >::setSrc().

00023                                                                                                                          :
00024   theRefitterAlgo( config ),
00025   theNavigator( navigator ),
00026   theDebugFlag( config.getUntrackedParameter<bool>( "debug", true ) )
00027 {
00028   TrackProducerBase< reco::Track >::setConf( config );
00029   TrackProducerBase< reco::Track >::setSrc( config.getParameter< edm::InputTag >( "src" ),
00030                                             config.getParameter< edm::InputTag >( "bsSrc" ) );
00031 }

KalmanAlignmentTrackRefitter::~KalmanAlignmentTrackRefitter ( void   ) 

Destructor.

Definition at line 34 of file KalmanAlignmentTrackRefitter.cc.

00034 {}


Member Function Documentation

void KalmanAlignmentTrackRefitter::debugTrackData ( const std::string  identifier,
const Trajectory traj,
const reco::Track track 
) [private]

Referenced by refitTracks().

virtual void KalmanAlignmentTrackRefitter::produce ( edm::Event ,
const edm::EventSetup  
) [inline, virtual]

Dummy implementation, due to inheritance from TrackProducerBase.

Implements TrackProducerBase< reco::Track >.

Definition at line 57 of file KalmanAlignmentTrackRefitter.h.

00057 {}

KalmanAlignmentTrackRefitter::TrajTrackPairCollection KalmanAlignmentTrackRefitter::refitSingleTracklet ( const TrackingGeometry geometry,
const MagneticField magneticField,
const TrajectoryFitter fitter,
const Propagator propagator,
const TransientTrackingRecHitBuilder recHitBuilder,
const reco::TransientTrack originalTrack,
RecHitContainer recHits,
const SortingDirection sortingDir,
bool  useExternalEstimate,
bool  reuseMomentumEstimate 
) [private]

Definition at line 169 of file KalmanAlignmentTrackRefitter.cc.

References anyDirection, asHepMatrix(), TransientTrackingRecHitBuilder::build(), reco::BeamSpot::dummy(), edm::OwnVector< T, P >::front(), reco::TransientTrack::impactPointState(), TrajectoryStateOnSurface::isValid(), it, TrajectoryStateOnSurface::localError(), TrajectoryStateOnSurface::localParameters(), LocalTrajectoryError::matrix(), AnalyticalPropagator::propagate(), Propagator::propagationDirection(), HLT_VtxMuL3::result, TrackProducerAlgorithm< T >::runWithCandidate(), edm::OwnVector< T, P >::size(), sortRecHits(), state, TrajectoryStateOnSurface::surface(), and theRefitterAlgo.

Referenced by refitTracks().

00179 {
00180 
00181   TrajTrackPairCollection result;
00182 
00183   if ( recHits.size() < 2 ) return result;
00184 
00185   sortRecHits( recHits, recHitBuilder, sortingDir );
00186 
00187   TransientTrackingRecHit::RecHitPointer firstHit = recHitBuilder->build( &recHits.front() );
00188 
00189   AnalyticalPropagator firstStatePropagator( magneticField, anyDirection );
00190   TrajectoryStateOnSurface firstState = firstStatePropagator.propagate( fullTrack.impactPointState(), firstHit->det()->surface() );
00191 
00192   if ( !firstState.isValid() ) return result;
00193 
00194   const double startErrorValue = 100;
00195   const unsigned int nTrajParam = 5;
00196 
00197   LocalTrajectoryError startError;
00198 
00199   if ( useExternalEstimate ) {
00200     startError = firstState.localError();
00201   } else {
00202     if ( reuseMomentumEstimate )
00203     {
00204       AlgebraicSymMatrix firstStateError( asHepMatrix( firstState.localError().matrix() ) );
00205       AlgebraicSymMatrix startErrorMatrix( nTrajParam, 0 );
00206       startErrorMatrix[0][0] = firstStateError[0][0];
00207       startErrorMatrix[1][1] = startErrorValue;//firstStateError[1][1];
00208       startErrorMatrix[2][2] = startErrorValue;//firstStateError[2][2];
00209       startErrorMatrix[3][3] = startErrorValue;
00210       startErrorMatrix[4][4] = startErrorValue;
00211       startError = LocalTrajectoryError( startErrorMatrix );
00212     } else {
00213       AlgebraicSymMatrix startErrorMatrix( startErrorValue*AlgebraicSymMatrix( nTrajParam, 1 ) );
00214       startError = LocalTrajectoryError( startErrorMatrix );
00215     }
00216 
00217   }
00218 
00219 //   // MOMENTUM ESTIMATE FOR COSMICS. P = 1.5 GeV
00220 //   LocalTrajectoryParameters firstStateParameters = firstState.localParameters();
00221 //   AlgebraicVector firstStateParamVec = asHepVector( firstStateParameters.mixedFormatVector() );
00222 //   firstStateParamVec[0] = 1./1.5;
00223 //   LocalTrajectoryParameters cosmicsStateParameters( firstStateParamVec, firstStateParameters.pzSign(), true );
00224 //   TrajectoryStateOnSurface tsos( cosmicsStateParameters, startError, firstState.surface(), magneticField );
00225 
00226   TrajectoryStateOnSurface tsos( firstState.localParameters(), startError, firstState.surface(), magneticField );
00227 
00228   TrajectoryStateTransform stateTransform;
00229   PTrajectoryStateOnDet state = *stateTransform.persistentState( tsos, firstHit->det()->geographicalId().rawId() );
00230 
00231   // Generate a trajectory seed.
00232   TrajectorySeed seed( state, recHits, propagator->propagationDirection() );
00233 
00234   // Generate track candidate.
00235   TrackCandidate candidate( recHits, seed, state );
00236 
00237   TrackCandidateCollection candidateCollection;
00238   candidateCollection.push_back( candidate );
00239 
00240   // Only dummy implementation of beam spot constraint.
00241   reco::BeamSpot dummyBeamSpot;
00242   dummyBeamSpot.dummy();
00243 
00244   AlgoProductCollection algoResult;
00245 
00246   theRefitterAlgo.runWithCandidate( geometry, magneticField, candidateCollection,
00247                                     fitter, propagator, recHitBuilder, dummyBeamSpot,
00248                                     algoResult );
00249 
00250   for ( AlgoProductCollection::iterator it = algoResult.begin(); it != algoResult.end(); ++it )
00251     result.push_back( make_pair( (*it).first, (*it).second.first ) );
00252 
00253   return result;
00254 }

KalmanAlignmentTrackRefitter::TrackletCollection KalmanAlignmentTrackRefitter::refitTracks ( const edm::EventSetup eventSetup,
const AlignmentSetupCollection algoSetups,
const ConstTrajTrackPairCollection tracks 
)

Definition at line 38 of file KalmanAlignmentTrackRefitter.cc.

References AlignableNavigator::alignableFromDetId(), anyDirection, debugTrackData(), edm::OwnVector< T, P >::empty(), TrackProducerBase< reco::Track >::getFromES(), edm::ESHandle< T >::product(), edm::OwnVector< T, P >::push_back(), refitSingleTracklet(), HLT_VtxMuL3::result, edm::OwnVector< T, P >::size(), theDebugFlag, theNavigator, and true.

Referenced by KalmanAlignmentAlgorithm::run().

00041 {
00042   // Retrieve what we need from the EventSetup
00043   edm::ESHandle< TrackerGeometry > aGeometry;
00044   edm::ESHandle< MagneticField > aMagneticField;
00045   edm::ESHandle< TrajectoryFitter > aTrajectoryFitter;
00046   edm::ESHandle< Propagator > aPropagator;
00047   edm::ESHandle< TransientTrackingRecHitBuilder > aRecHitBuilder;
00048 
00049   getFromES( setup, aGeometry, aMagneticField, aTrajectoryFitter, aPropagator, aRecHitBuilder );
00050 
00051   TrackletCollection result;
00052   TrackCollection fullTracks;
00053 
00054   ConstTrajTrackPairCollection refittedFullTracks;
00055   ConstTrajTrackPairCollection::const_iterator itTrack;
00056 
00057   for( itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack )
00058   {
00059     TransientTrack fullTrack( *(*itTrack).second, aMagneticField.product() );
00060 
00061     AlignmentSetupCollection::const_iterator itSetup;
00062     for ( itSetup = algoSetups.begin(); itSetup != algoSetups.end(); ++itSetup )
00063     {
00064       RecHitContainer trackingRecHits;
00065       RecHitContainer externalTrackingRecHits;
00066 
00067       // Extract collection with TrackingRecHits
00068       Trajectory::ConstRecHitContainer hits = (*itTrack).first->recHits();
00069       Trajectory::ConstRecHitContainer::iterator itHits;
00070 
00071       for ( itHits = hits.begin(); itHits != hits.end(); ++itHits )
00072       {
00073         if ( !(*itHits)->isValid() ) continue;
00074 
00075         try
00076         {
00077           //if ( !theNavigator->alignableFromDetId( (*itHits)->geographicalId() )->alignmentParameters() ) continue;
00078           theNavigator->alignableFromDetId( (*itHits)->geographicalId() );        
00079         } catch(...) { continue; }
00080 
00081         if ( (*itSetup)->useForTracking( *itHits ) )
00082         {
00083           trackingRecHits.push_back( (*itHits)->hit()->clone() );
00084         }
00085         else if ( (*itSetup)->useForExternalTracking( *itHits ) )
00086         {
00087           externalTrackingRecHits.push_back( (*itHits)->hit()->clone() );
00088         }
00089       }
00090 
00091       //edm::LogInfo( "KalmanAlignmentTrackRefitter" ) << "Hits for tracking/external: " << trackingRecHits.size() << "/" << externalTrackingRecHits.size();
00092 
00093       if ( trackingRecHits.empty() ) continue;
00094 
00095       if ( externalTrackingRecHits.empty() )
00096       {
00097         if ( ( (*itSetup)->getExternalTrackingSubDetIds().size() == 0 ) && // O.K., no external hits expected,
00098              ( trackingRecHits.size() >= (*itSetup)->minTrackingHits() ) )
00099         {
00100           TrajTrackPairCollection refitted = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
00101                                                                   (*itSetup)->fitter(), (*itSetup)->propagator(), 
00102                                                                   aRecHitBuilder.product(), fullTrack,
00103                                                                   trackingRecHits, (*itSetup)->sortingDirection(),
00104                                                                   false, true );
00105 
00106           if ( refitted.empty() ) continue; // The refitting did not work ... Try next!
00107 
00108           if ( theDebugFlag )
00109           {
00110             debugTrackData( (*itSetup)->id(), refitted.front().first, refitted.front().second );
00111             debugTrackData( "OrigFullTrack", (*itTrack).first, (*itTrack).second );
00112           }
00113 
00114 
00115           TrackletPtr trackletPtr( new KalmanAlignmentTracklet( refitted.front(), *itSetup ) );
00116           result.push_back( trackletPtr );
00117         }
00118         else { continue; } // Expected external hits but found none or not enough hits.
00119       }
00120       else if ( ( trackingRecHits.size() >= (*itSetup)->minTrackingHits() ) &&
00121                 ( externalTrackingRecHits.size() >= (*itSetup)->minExternalHits() ) ) 
00122       {
00123         // Create an instance of KalmanAlignmentTracklet with an external prediction.
00124 
00125         TrajTrackPairCollection external = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
00126                                                                 (*itSetup)->externalFitter(), (*itSetup)->externalPropagator(),
00127                                                                 aRecHitBuilder.product(), fullTrack,
00128                                                                 externalTrackingRecHits, (*itSetup)->externalSortingDirection(),
00129                                                                 false, true );
00130         if ( external.empty() ) { continue; }
00131 
00132         TransientTrack externalTrack( *external.front().second, aMagneticField.product() );
00133 
00134         TrajTrackPairCollection refitted = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
00135                                                                 (*itSetup)->fitter(), (*itSetup)->propagator(),
00136                                                                 aRecHitBuilder.product(), externalTrack,
00137                                                                 trackingRecHits, (*itSetup)->sortingDirection(),
00138                                                                 false, true );
00139         if ( refitted.empty() ) continue;
00140 
00141         //const Surface& surface = refitted.front().first->firstMeasurement().updatedState().surface();
00142         const Surface& surface = refitted.front().first->lastMeasurement().updatedState().surface();
00143         TrajectoryStateOnSurface externalTsos = externalTrack.impactPointState();
00144         AnalyticalPropagator externalPredictionPropagator( aMagneticField.product(), anyDirection );
00145         TrajectoryStateOnSurface externalPrediction = externalPredictionPropagator.propagate( externalTsos, surface );
00146         if ( !externalPrediction.isValid() ) continue;
00147 
00148         if ( theDebugFlag )
00149         {
00150           debugTrackData( string("External") + (*itSetup)->id(), external.front().first, external.front().second );
00151           debugTrackData( (*itSetup)->id(), refitted.front().first, refitted.front().second );
00152           debugTrackData( "OrigFullTrack", (*itTrack).first, (*itTrack).second );
00153         }
00154 
00155         TrackletPtr trackletPtr( new KalmanAlignmentTracklet( refitted.front(), externalPrediction, *itSetup ) );
00156         result.push_back( trackletPtr );
00157 
00158         delete external.front().first;
00159         delete external.front().second;
00160       }
00161     }
00162   }
00163 
00164   return result;
00165 }

void KalmanAlignmentTrackRefitter::sortRecHits ( RecHitContainer hits,
const TransientTrackingRecHitBuilder builder,
const SortingDirection sortingDir 
) const [private]

Definition at line 257 of file KalmanAlignmentTrackRefitter.cc.

References edm::OwnVector< T, P >::back(), edm::OwnVector< T, P >::begin(), TransientTrackingRecHitBuilder::build(), edm::OwnVector< T, P >::end(), edm::OwnVector< T, P >::front(), insideOut, muonGeometry::mag(), edm::OwnVector< T, P >::push_back(), edm::OwnVector< T, P >::size(), KalmanAlignmentSetup::sortDownsideUp, KalmanAlignmentSetup::sortInsideOut, KalmanAlignmentSetup::sortOutsideIn, KalmanAlignmentSetup::sortUpsideDown, edm::OwnVector< T, P >::swap(), tmp, and y.

Referenced by refitSingleTracklet().

00260 {
00261   // Don't start sorting if there is only 1 or even 0 elements.
00262   if ( hits.size() < 2 ) return;
00263 
00264   TransientTrackingRecHit::RecHitPointer firstHit = builder->build( &hits.front() );
00265   double firstRadius = firstHit->det()->surface().toGlobal( firstHit->localPosition() ).mag();
00266   double firstY = firstHit->det()->surface().toGlobal( firstHit->localPosition() ).y();
00267 
00268   TransientTrackingRecHit::RecHitPointer lastHit = builder->build( &hits.back() );
00269   double lastRadius = lastHit->det()->surface().toGlobal( lastHit->localPosition() ).mag();
00270   double lastY = lastHit->det()->surface().toGlobal( lastHit->localPosition() ).y();
00271 
00272   bool insideOut = firstRadius < lastRadius;
00273   bool upsideDown = lastY < firstY;
00274 
00275   if ( ( insideOut && ( sortingDir == KalmanAlignmentSetup::sortInsideOut ) ) ||
00276        ( !insideOut && ( sortingDir == KalmanAlignmentSetup::sortOutsideIn ) ) ||
00277        ( upsideDown && ( sortingDir == KalmanAlignmentSetup::sortUpsideDown ) ) ||
00278        ( !upsideDown && ( sortingDir == KalmanAlignmentSetup::sortDownsideUp ) ) ) return;
00279 
00280   // Fill temporary container with reversed hits.
00281   RecHitContainer tmp;
00282   RecHitContainer::iterator itHit = hits.end();
00283   do { --itHit; tmp.push_back( ( *itHit ).clone() ); } while ( itHit != hits.begin() );
00284 
00285   // Swap the content of the temporary and the input container.
00286   hits.swap( tmp );
00287 
00288   return;
00289 }


Member Data Documentation

bool KalmanAlignmentTrackRefitter::theDebugFlag [private]

Definition at line 80 of file KalmanAlignmentTrackRefitter.h.

Referenced by refitTracks().

AlignableNavigator* KalmanAlignmentTrackRefitter::theNavigator [private]

Definition at line 79 of file KalmanAlignmentTrackRefitter.h.

Referenced by refitTracks().

TrackProducerAlgorithm<reco::Track> KalmanAlignmentTrackRefitter::theRefitterAlgo [private]

Definition at line 78 of file KalmanAlignmentTrackRefitter.h.

Referenced by refitSingleTracklet().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:26:15 2009 for CMSSW by  doxygen 1.5.4