CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/Alignment/KalmanAlignmentAlgorithm/src/KalmanAlignmentTrackRefitter.cc

Go to the documentation of this file.
00001 
00002 #include "Alignment/KalmanAlignmentAlgorithm/interface/KalmanAlignmentTrackRefitter.h"
00003 
00004 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00005 
00006 #include "RecoTracker/TransientTrackingRecHit/interface/TRecHit1DMomConstraint.h"
00007 
00008 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00009 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00010 
00011 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00012 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00013 
00014 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00015 
00016 #include "Alignment/KalmanAlignmentAlgorithm/interface/KalmanAlignmentDataCollector.h"
00017 
00018 #include "CLHEP/GenericFunctions/CumulativeChiSquare.hh"
00019 
00020 #include <iostream>
00021 
00022 using namespace std;
00023 using namespace reco;
00024 using namespace Genfun;
00025 
00026 
00027 KalmanAlignmentTrackRefitter::KalmanAlignmentTrackRefitter( const edm::ParameterSet& config, AlignableNavigator* navigator ) :
00028   theRefitterAlgo( config ),
00029   theNavigator( navigator ),
00030   theDebugFlag( config.getUntrackedParameter<bool>( "debug", true ) )
00031 {
00032   TrackProducerBase< reco::Track >::setConf( config );
00033   TrackProducerBase< reco::Track >::setSrc( config.getParameter< edm::InputTag >( "src" ),
00034                                             config.getParameter< edm::InputTag >( "bsSrc" ) );
00035 }
00036 
00037 
00038 KalmanAlignmentTrackRefitter::~KalmanAlignmentTrackRefitter( void ) {}
00039 
00040 
00041 KalmanAlignmentTrackRefitter::TrackletCollection
00042 KalmanAlignmentTrackRefitter::refitTracks( const edm::EventSetup& setup,
00043                                            const AlignmentSetupCollection& algoSetups,
00044                                            const ConstTrajTrackPairCollection& tracks,
00045                                            const reco::BeamSpot* beamSpot )
00046 {
00047   // Retrieve what we need from the EventSetup
00048   edm::ESHandle< TrackerGeometry > aGeometry;
00049   edm::ESHandle< MagneticField > aMagneticField;
00050   edm::ESHandle< TrajectoryFitter > aTrajectoryFitter;
00051   edm::ESHandle< Propagator > aPropagator;
00052   edm::ESHandle<MeasurementTracker> theMeasTk;
00053   edm::ESHandle< TransientTrackingRecHitBuilder > aRecHitBuilder;
00054 
00055   getFromES( setup, aGeometry, aMagneticField, aTrajectoryFitter, aPropagator, theMeasTk, aRecHitBuilder );
00056 
00057   TrackletCollection result;
00058   TrackCollection fullTracks;
00059 
00060   ConstTrajTrackPairCollection refittedFullTracks;
00061   ConstTrajTrackPairCollection::const_iterator itTrack;
00062 
00063   for( itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack )
00064   {
00065     TransientTrack fullTrack( *(*itTrack).second, aMagneticField.product() );
00066 
00067     AlignmentSetupCollection::const_iterator itSetup;
00068     for ( itSetup = algoSetups.begin(); itSetup != algoSetups.end(); ++itSetup )
00069     {
00070       RecHitContainer trackingRecHits;
00071       RecHitContainer externalTrackingRecHits;
00072 
00074       RecHitContainer zPlusRecHits;
00075       RecHitContainer zMinusRecHits;
00076 
00077       // Extract collection with TrackingRecHits
00078       Trajectory::ConstRecHitContainer hits = (*itTrack).first->recHits();
00079       Trajectory::ConstRecHitContainer::iterator itHits;
00080 
00081       for ( itHits = hits.begin(); itHits != hits.end(); ++itHits )
00082       {
00083         if ( !(*itHits)->isValid() ) continue;
00084 
00085         try
00086         {
00087           //if ( !theNavigator->alignableFromDetId( (*itHits)->geographicalId() )->alignmentParameters() ) continue;
00088           theNavigator->alignableFromDetId( (*itHits)->geographicalId() );        
00089         } catch(...) { continue; }
00090 
00091         if ( (*itSetup)->useForTracking( *itHits ) )
00092         {
00093           trackingRecHits.push_back( (*itHits)->hit()->clone() );
00094 
00095           ( (*itHits)->det()->position().z() > 0. ) ?
00096             zPlusRecHits.push_back( (*itHits)->hit()->clone() ) :
00097             zMinusRecHits.push_back( (*itHits)->hit()->clone() );
00098 
00101         }
00102         else if ( (*itSetup)->useForExternalTracking( *itHits ) )
00103         {
00104           externalTrackingRecHits.push_back( (*itHits)->hit()->clone() );
00105         }
00106       }
00107 
00108       //edm::LogInfo( "KalmanAlignmentTrackRefitter" ) << "Hits for tracking/external: " << trackingRecHits.size() << "/" << externalTrackingRecHits.size();
00109 
00110       //if ( !zPlusRecHits.size() || !zMinusRecHits.size() ) continue;
00112 
00113       if ( trackingRecHits.empty() ) continue;
00114 
00115       if ( externalTrackingRecHits.empty() )
00116       {
00117         if ( ( (*itSetup)->getExternalTrackingSubDetIds().size() == 0 ) && // O.K., no external hits expected,
00118              ( trackingRecHits.size() >= (*itSetup)->minTrackingHits() ) )
00119         {
00120           TrajTrackPairCollection refitted = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
00121                                                                   (*itSetup)->fitter(), (*itSetup)->propagator(), 
00122                                                                   aRecHitBuilder.product(), fullTrack,
00123                                                                   trackingRecHits, beamSpot,
00124                                                                   (*itSetup)->sortingDirection(), false, true );
00125 
00126           // The refitting did not work ... Try next!
00127           if ( refitted.empty() ) continue;
00128           if ( rejectTrack( refitted.front().second ) ) continue;
00129 
00130           if ( theDebugFlag )
00131           {
00132             debugTrackData( (*itSetup)->id(), refitted.front().first, refitted.front().second, beamSpot );
00133             debugTrackData( "OrigFullTrack", (*itTrack).first, (*itTrack).second, beamSpot );
00134           }
00135 
00136 
00137           TrackletPtr trackletPtr( new KalmanAlignmentTracklet( refitted.front(), *itSetup ) );
00138           result.push_back( trackletPtr );
00139         }
00140         else { continue; } // Expected external hits but found none or not enough hits.
00141       }
00142       else if ( ( trackingRecHits.size() >= (*itSetup)->minTrackingHits() ) &&
00143                 ( externalTrackingRecHits.size() >= (*itSetup)->minExternalHits() ) ) 
00144       {
00145         // Create an instance of KalmanAlignmentTracklet with an external prediction.
00146 
00147         TrajTrackPairCollection external = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
00148                                                                 (*itSetup)->externalFitter(), (*itSetup)->externalPropagator(),
00149                                                                 aRecHitBuilder.product(), fullTrack,
00150                                                                 externalTrackingRecHits, beamSpot,
00151                                                                 (*itSetup)->externalSortingDirection(),
00152                                                                 false, true );
00153         //if ( external.empty() || rejectTrack( external.front().second ) ) { continue; }
00154         if ( external.empty() ) { continue; }
00155 
00156         TransientTrack externalTrack( *external.front().second, aMagneticField.product() );
00157 
00158         TrajTrackPairCollection refitted = refitSingleTracklet( aGeometry.product(), aMagneticField.product(),
00159                                                                 (*itSetup)->fitter(), (*itSetup)->propagator(),
00160                                                                 aRecHitBuilder.product(), externalTrack,
00161                                                                 trackingRecHits, beamSpot,
00162                                                                 (*itSetup)->sortingDirection(),
00163                                                                 false, true, (*itSetup)->id() );
00164 
00165         if ( refitted.empty() ) { continue; }
00166         if ( rejectTrack( refitted.front().second ) ) continue;
00167 
00168         //const Surface& surface = refitted.front().first->firstMeasurement().updatedState().surface();
00169         const Surface& surface = refitted.front().first->lastMeasurement().updatedState().surface();
00170         TrajectoryStateOnSurface externalTsos = externalTrack.impactPointState();
00171         AnalyticalPropagator externalPredictionPropagator( aMagneticField.product(), anyDirection );
00172         TrajectoryStateOnSurface externalPrediction = externalPredictionPropagator.propagate( externalTsos, surface );
00173         if ( !externalPrediction.isValid() ) continue;
00174 
00175         if ( theDebugFlag )
00176         {
00177           debugTrackData( string("External") + (*itSetup)->id(), external.front().first, external.front().second, beamSpot );
00178           debugTrackData( (*itSetup)->id(), refitted.front().first, refitted.front().second, beamSpot );
00179           debugTrackData( "OrigFullTrack", (*itTrack).first, (*itTrack).second, beamSpot );
00180         }
00181 
00182         TrackletPtr trackletPtr( new KalmanAlignmentTracklet( refitted.front(), externalPrediction, *itSetup ) );
00183         result.push_back( trackletPtr );
00184 
00185         delete external.front().first;
00186         delete external.front().second;
00187       }
00188     }
00189   }
00190 
00191   return result;
00192 }
00193 
00194 
00195 KalmanAlignmentTrackRefitter::TrajTrackPairCollection
00196 KalmanAlignmentTrackRefitter::refitSingleTracklet( const TrackingGeometry* geometry,
00197                                                    const MagneticField* magneticField,
00198                                                    const TrajectoryFitter* fitter,
00199                                                    const Propagator* propagator,
00200                                                    const TransientTrackingRecHitBuilder* recHitBuilder,
00201                                                    const TransientTrack& fullTrack,
00202                                                    RecHitContainer& recHits,
00203                                                    const reco::BeamSpot* beamSpot,
00204                                                    const SortingDirection& sortingDir,
00205                                                    bool useExternalEstimate,
00206                                                    bool reuseMomentumEstimate,
00207                                                    const string identifier )
00208 {
00209 
00210   TrajTrackPairCollection result;
00211 
00212   if ( recHits.size() < 2 ) return result;
00213 
00214   sortRecHits( recHits, recHitBuilder, sortingDir );
00215 
00216   TransientTrackingRecHit::RecHitContainer hits;
00217   RecHitContainer::iterator itRecHit;
00218   for ( itRecHit = recHits.begin(); itRecHit != recHits.end(); ++itRecHit )
00219     hits.push_back( recHitBuilder->build( &(*itRecHit) ) );
00220 
00221   TransientTrackingRecHit::ConstRecHitPointer firstHit = hits.front();
00222 
00223   AnalyticalPropagator firstStatePropagator( magneticField, anyDirection );
00224   TrajectoryStateOnSurface firstState = firstStatePropagator.propagate( fullTrack.impactPointState(), firstHit->det()->surface() );
00225 
00226   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_IPPt"), 
00227                                                1e-2*fullTrack.impactPointState().globalParameters().momentum().perp() );
00228 
00229   if ( !firstState.isValid() ) return result;
00230 
00231 //   LocalTrajectoryError startError;
00232 
00233 //   const double startErrorValue = 100;
00234 //   const unsigned int nTrajParam = 5;
00235 
00236 //   if ( useExternalEstimate ) {
00237 //     startError = firstState.localError();
00238 //   } else {
00239 //     if ( reuseMomentumEstimate )
00240 //     {
00241 //       AlgebraicSymMatrix firstStateError( asHepMatrix( firstState.localError().matrix() ) );
00242 //       AlgebraicSymMatrix startErrorMatrix( nTrajParam, 0 );
00243 //       startErrorMatrix[0][0] = 1e-10;
00244 //       //startErrorMatrix[0][0] = firstStateError[0][0];
00245 //       startErrorMatrix[1][1] = startErrorValue;//firstStateError[1][1];
00246 //       startErrorMatrix[2][2] = startErrorValue;//firstStateError[2][2];
00247 //       startErrorMatrix[3][3] = startErrorValue;
00248 //       startErrorMatrix[4][4] = startErrorValue;
00249 //       startError = LocalTrajectoryError( startErrorMatrix );
00250 //     } else {
00251 //       AlgebraicSymMatrix startErrorMatrix( startErrorValue*AlgebraicSymMatrix( nTrajParam, 1 ) );
00252 //       startError = LocalTrajectoryError( startErrorMatrix );
00253 //     }
00254 
00255 //   }
00256 
00257 //   // MOMENTUM ESTIMATE FOR COSMICS. P = 1.5 GeV
00258 //   LocalTrajectoryParameters firstStateParameters = firstState.localParameters();
00259 //   AlgebraicVector firstStateParamVec = asHepVector( firstStateParameters.mixedFormatVector() );
00260 //   firstStateParamVec[0] = 1./1.5;
00261 //   LocalTrajectoryParameters cosmicsStateParameters( firstStateParamVec, firstStateParameters.pzSign(), true );
00262 //   TrajectoryStateOnSurface tsos( cosmicsStateParameters, startError, firstState.surface(), magneticField );
00263 
00264   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_FSPt"),
00265                                                1e-2*firstState.globalParameters().momentum().perp() );
00266 
00267   firstState.rescaleError( 100 );
00268   TrajectoryStateOnSurface tsos( firstState.localParameters(), firstState.localError(),
00269                                  firstState.surface(), magneticField );
00270 
00271   // Generate a trajectory seed.
00272   TrajectorySeed seed( PTrajectoryStateOnDet(), recHits, propagator->propagationDirection() );
00273 
00274   // Generate track candidate.
00275   
00276   PTrajectoryStateOnDet state = trajectoryStateTransform::persistentState( tsos, firstHit->det()->geographicalId().rawId() );
00277   TrackCandidate candidate( recHits, seed, state );
00278 
00279   AlgoProductCollection algoResult;
00280 
00281   int charge = static_cast<int>( tsos.charge() );
00282   double momentum = firstState.localParameters().momentum().mag();
00283   TransientTrackingRecHit::RecHitPointer testhit =
00284     TRecHit1DMomConstraint::build( charge, momentum, 1e-10, &tsos.surface() );
00285 
00286   //no insert in OwnVector...
00287   TransientTrackingRecHit::RecHitContainer tmpHits;
00288   tmpHits.push_back(testhit);
00289   for (TransientTrackingRecHit::RecHitContainer::const_iterator i=hits.begin(); i!=hits.end(); i++){
00290     tmpHits.push_back(*i);
00291   }
00292   hits.swap(tmpHits);
00293 
00294   theRefitterAlgo.buildTrack( fitter, propagator, algoResult, hits, tsos, seed, 0, *beamSpot, candidate.seedRef());
00295 
00296   for ( AlgoProductCollection::iterator it = algoResult.begin(); it != algoResult.end(); ++it )
00297     result.push_back( make_pair( (*it).first, (*it).second.first ) );
00298 
00299   return result;
00300 }
00301 
00302 
00303 void KalmanAlignmentTrackRefitter::sortRecHits( RecHitContainer& hits,
00304                                                 const TransientTrackingRecHitBuilder* builder,
00305                                                 const SortingDirection& sortingDir ) const
00306 {
00307   // Don't start sorting if there is only 1 or even 0 elements.
00308   if ( hits.size() < 2 ) return;
00309 
00310   TransientTrackingRecHit::RecHitPointer firstHit = builder->build( &hits.front() );
00311   double firstRadius = firstHit->det()->surface().toGlobal( firstHit->localPosition() ).mag();
00312   double firstY = firstHit->det()->surface().toGlobal( firstHit->localPosition() ).y();
00313 
00314   TransientTrackingRecHit::RecHitPointer lastHit = builder->build( &hits.back() );
00315   double lastRadius = lastHit->det()->surface().toGlobal( lastHit->localPosition() ).mag();
00316   double lastY = lastHit->det()->surface().toGlobal( lastHit->localPosition() ).y();
00317 
00318   bool insideOut = firstRadius < lastRadius;
00319   bool upsideDown = lastY < firstY;
00320 
00321   if ( ( insideOut && ( sortingDir == KalmanAlignmentSetup::sortInsideOut ) ) ||
00322        ( !insideOut && ( sortingDir == KalmanAlignmentSetup::sortOutsideIn ) ) ||
00323        ( upsideDown && ( sortingDir == KalmanAlignmentSetup::sortUpsideDown ) ) ||
00324        ( !upsideDown && ( sortingDir == KalmanAlignmentSetup::sortDownsideUp ) ) ) return;
00325 
00326   // Fill temporary container with reversed hits.
00327   RecHitContainer tmp;
00328   RecHitContainer::iterator itHit = hits.end();
00329   do { --itHit; tmp.push_back( ( *itHit ).clone() ); } while ( itHit != hits.begin() );
00330 
00331   // Swap the content of the temporary and the input container.
00332   hits.swap( tmp );
00333 
00334   return;
00335 }
00336 
00337 
00338 bool KalmanAlignmentTrackRefitter::rejectTrack( const Track* track ) const
00339 {
00340   double trackChi2 = track->chi2();
00341   unsigned int ndof = static_cast<unsigned int>( track->ndof() );
00342   if ( trackChi2 <= 0. || ndof <= 0 ) return false;
00343 
00344   //FIXME: should be configurable (via KalmanAlignmentSetup)
00345   double minChi2Prob = 0;//1e-6;
00346   double maxChi2Prob = 1.0;
00347 
00348   GENFUNCTION cumulativeChi2 = Genfun::CumulativeChiSquare( ndof );
00349   double chi2Prob = 1. - cumulativeChi2( trackChi2 );
00350   return ( chi2Prob < minChi2Prob ) || ( chi2Prob > maxChi2Prob ); 
00351 }
00352 
00353 
00354 void KalmanAlignmentTrackRefitter::debugTrackData( const string identifier,
00355                                                    const Trajectory* traj,
00356                                                    const Track* track,
00357                                                    const reco::BeamSpot* bs )
00358 {
00359   unsigned int ndof = static_cast<unsigned int>( track->ndof() );
00360   double trackChi2 = track->chi2();
00361   if ( ( trackChi2 > 0. ) && ( ndof > 0 ) )
00362   {
00363     GENFUNCTION cumulativeChi2 = Genfun::CumulativeChiSquare( ndof );
00364     KalmanAlignmentDataCollector::fillHistogram( identifier + string("_CumChi2"), 1. - cumulativeChi2( trackChi2 ) );
00365   } else if ( ndof == 0 ) {
00366     KalmanAlignmentDataCollector::fillHistogram( identifier + string("_CumChi2"), -1. );
00367   } else {
00368     KalmanAlignmentDataCollector::fillHistogram( identifier + string("_CumChi2"), -2. );
00369   }
00370 
00371   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_NHits"), traj->foundHits() );
00372   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_Pt"), 1e-2*track->pt() );
00373   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_Eta"), track->eta() );
00374   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_Phi"), track->phi() );
00375   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_NormChi2"), track->normalizedChi2() );
00376   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_DZ"), track->dz() );
00377 
00378   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_DXY_BS"), fabs( track->dxy( bs->position() ) ) );
00379   KalmanAlignmentDataCollector::fillHistogram( identifier + string("_DXY"), fabs( track->dxy() ) );
00380   //KalmanAlignmentDataCollector::fillHistogram( identifier + string("_D0"), fabs( track->d0() ) );
00381 }