CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/RecoTracker/CkfPattern/src/CkfTrajectoryBuilder.cc

Go to the documentation of this file.
00001 #include "RecoTracker/CkfPattern/interface/CkfTrajectoryBuilder.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00006 
00007 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00008 #include "TrackingTools/PatternTools/interface/TrajectoryStateUpdator.h"
00009 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimatorBase.h"
00010 
00011 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00012 #include "TrackingTools/PatternTools/interface/TrajMeasLessEstim.h"
00013 #include "TrackingTools/TrajectoryState/interface/BasicSingleTrajectoryState.h"
00014 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00015 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00016 
00017 
00018 #include "RecoTracker/CkfPattern/src/RecHitIsInvalid.h"
00019 #include "RecoTracker/CkfPattern/interface/TrajCandLess.h"
00020 
00021 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00022 
00023 #include "RecoTracker/CkfPattern/interface/IntermediateTrajectoryCleaner.h"
00024 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00025 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00026 
00027 #include "FWCore/Services/interface/UpdaterService.h"
00028 
00029 using namespace std;
00030 
00031 
00032 CkfTrajectoryBuilder::
00033   CkfTrajectoryBuilder(const edm::ParameterSet&              conf,
00034                        const TrajectoryStateUpdator*         updator,
00035                        const Propagator*                     propagatorAlong,
00036                        const Propagator*                     propagatorOpposite,
00037                        const Chi2MeasurementEstimatorBase*   estimator,
00038                        const TransientTrackingRecHitBuilder* recHitBuilder,
00039                        const MeasurementTracker*             measurementTracker,
00040                        const TrajectoryFilter*               filter):
00041 
00042     BaseCkfTrajectoryBuilder(conf,
00043                              updator, propagatorAlong,propagatorOpposite,
00044                              estimator, recHitBuilder, measurementTracker,filter)
00045 {
00046   theMaxCand              = conf.getParameter<int>("maxCand");
00047   theLostHitPenalty       = conf.getParameter<double>("lostHitPenalty");
00048   theIntermediateCleaning = conf.getParameter<bool>("intermediateCleaning");
00049   theAlwaysUseInvalidHits = conf.getParameter<bool>("alwaysUseInvalidHits");
00050   /*
00051     theSharedSeedCheck = conf.getParameter<bool>("SharedSeedCheck");
00052     std::stringstream ss;
00053     ss<<"CkfTrajectoryBuilder_"<<conf.getParameter<std::string>("ComponentName")<<"_"<<this;
00054     theUniqueName = ss.str();
00055     LogDebug("CkfPattern")<<"my unique name is: "<<theUniqueName;
00056   */
00057 }
00058 
00059 /*
00060   void CkfTrajectoryBuilder::setEvent(const edm::Event& event) const
00061   {
00062   theMeasurementTracker->update(event);
00063   }
00064 */
00065 
00066 CkfTrajectoryBuilder::TrajectoryContainer 
00067 CkfTrajectoryBuilder::trajectories(const TrajectorySeed& seed) const
00068 {  
00069   TrajectoryContainer result;
00070   result.reserve(5);
00071   trajectories(seed, result);
00072   return result;
00073 }
00074 
00075 /*
00076   void CkfTrajectoryBuilder::rememberSeedAndTrajectories(const TrajectorySeed& seed,
00077   CkfTrajectoryBuilder::TrajectoryContainer &result) const
00078   {
00079   
00080   //result ----> theCachedTrajectories
00081   //every first iteration on event. forget about everything that happened before
00082   if (edm::Service<UpdaterService>()->checkOnce(theUniqueName)) 
00083   theCachedTrajectories.clear();
00084   
00085   if (seed.nHits()==0) return;
00086   
00087   //then remember those trajectories
00088   for (TrajectoryContainer::iterator traj=result.begin();
00089   traj!=result.end(); ++traj) {
00090   theCachedTrajectories.insert(std::make_pair(seed.recHits().first->geographicalId(),*traj));
00091   }  
00092   }
00093   
00094   bool CkfTrajectoryBuilder::sharedSeed(const TrajectorySeed& s1,const TrajectorySeed& s2) const{
00095   //quit right away on nH=0
00096   if (s1.nHits()==0 || s2.nHits()==0) return false;
00097   //quit right away if not the same number of hits
00098   if (s1.nHits()!=s2.nHits()) return false;
00099   TrajectorySeed::range r1=s1.recHits();
00100   TrajectorySeed::range r2=s2.recHits();
00101   TrajectorySeed::const_iterator i1,i2;
00102   TrajectorySeed::const_iterator & i1_e=r1.second,&i2_e=r2.second;
00103   TrajectorySeed::const_iterator & i1_b=r1.first,&i2_b=r2.first;
00104   //quit right away if first detId does not match. front exist because of ==0 ->quit test
00105   if(i1_b->geographicalId() != i2_b->geographicalId()) return false;
00106   //then check hit by hit if they are the same
00107   for (i1=i1_b,i2=i2_b;i1!=i1_e,i2!=i2_e;++i1,++i2){
00108   if (!i1->sharesInput(&(*i2),TrackingRecHit::all)) return false;
00109   }
00110   return true;
00111   }
00112   bool CkfTrajectoryBuilder::seedAlreadyUsed(const TrajectorySeed& seed,
00113   CkfTrajectoryBuilder::TempTrajectoryContainer &candidates) const
00114   {
00115   //theCachedTrajectories ---> candidates
00116   if (seed.nHits()==0) return false;
00117   bool answer=false;
00118   pair<SharedTrajectory::const_iterator, SharedTrajectory::const_iterator> range = 
00119   theCachedTrajectories.equal_range(seed.recHits().first->geographicalId());
00120   SharedTrajectory::const_iterator trajP;
00121   for (trajP = range.first; trajP!=range.second;++trajP){
00122   //check whether seeds are identical     
00123   if (sharedSeed(trajP->second.seed(),seed)){
00124   candidates.push_back(trajP->second);
00125   answer=true;
00126   }//already existing trajectory shares the seed.   
00127   }//loop already made trajectories      
00128   
00129   return answer;
00130   }
00131 */
00132 
00133 void
00134 CkfTrajectoryBuilder::trajectories(const TrajectorySeed& seed, CkfTrajectoryBuilder::TrajectoryContainer &result) const
00135 {  
00136   // analyseSeed( seed);
00137   /*
00138     if (theSharedSeedCheck){
00139     TempTrajectoryContainer candidates;
00140     if (seedAlreadyUsed(seed,candidates))
00141     {
00142     //start with those candidates already made before
00143     limitedCandidates(candidates,result);
00144     //and quit
00145     return;
00146     }
00147     }
00148   */
00149 
00150   TempTrajectory startingTraj = createStartingTrajectory( seed );
00151 
00154   limitedCandidates( startingTraj, result);
00155 
00156   /*
00157     //and remember what you just did
00158     if (theSharedSeedCheck)  rememberSeedAndTrajectories(seed,result);
00159   */
00160 
00161   // analyseResult(result);
00162 }
00163 
00164 void CkfTrajectoryBuilder::
00165 limitedCandidates( TempTrajectory& startingTraj,
00166                    TrajectoryContainer& result) const
00167 {
00168   TempTrajectoryContainer candidates;
00169   candidates.push_back( startingTraj);
00170   limitedCandidates(candidates,result);
00171 }
00172 
00173 void CkfTrajectoryBuilder::
00174 limitedCandidates( TempTrajectoryContainer &candidates,
00175                    TrajectoryContainer& result) const
00176 {
00177   unsigned int nIter=1;
00178   //  TempTrajectoryContainer candidates; // = TrajectoryContainer();
00179   TempTrajectoryContainer newCand; // = TrajectoryContainer();
00180   //  candidates.push_back( startingTraj);
00181 
00182   while ( !candidates.empty()) {
00183 
00184     newCand.clear();
00185     for (TempTrajectoryContainer::iterator traj=candidates.begin();
00186          traj!=candidates.end(); traj++) {
00187       std::vector<TM> meas;
00188       findCompatibleMeasurements(*traj, meas);
00189 
00190       // --- method for debugging
00191       if(!analyzeMeasurementsDebugger(*traj,meas,
00192                                       theMeasurementTracker,
00193                                       theForwardPropagator,theEstimator,
00194                                       theTTRHBuilder)) return;
00195       // ---
00196 
00197       if ( meas.empty()) {
00198         if ( qualityFilter( *traj)) addToResult( *traj, result);
00199       }
00200       else {
00201         std::vector<TM>::const_iterator last;
00202         if ( theAlwaysUseInvalidHits) last = meas.end();
00203         else {
00204           if (meas.front().recHit()->isValid()) {
00205             last = find_if( meas.begin(), meas.end(), RecHitIsInvalid());
00206           }
00207           else last = meas.end();
00208         }
00209 
00210         for( std::vector<TM>::const_iterator itm = meas.begin(); 
00211              itm != last; itm++) {
00212           TempTrajectory newTraj = *traj;
00213           updateTrajectory( newTraj, *itm);
00214 
00215           if ( toBeContinued(newTraj)) {
00216             newCand.push_back(newTraj);
00217           }
00218           else {
00219             if ( qualityFilter(newTraj)) addToResult( newTraj, result);
00221           }
00222         }
00223       }
00224 
00225       if ((int)newCand.size() > theMaxCand) {
00226         sort( newCand.begin(), newCand.end(), TrajCandLess<TempTrajectory>(theLostHitPenalty));
00227         newCand.erase( newCand.begin()+theMaxCand, newCand.end());
00228       }
00229     }
00230 
00231     if (theIntermediateCleaning) IntermediateTrajectoryCleaner::clean(newCand);
00232 
00233     candidates.swap(newCand);
00234     
00235     LogDebug("CkfPattern") <<result.size()<<" candidates after "<<nIter++<<" CKF iteration: \n"
00236                            <<PrintoutHelper::dumpCandidates(result)
00237                            <<"\n "<<candidates.size()<<" running candidates are: \n"
00238                            <<PrintoutHelper::dumpCandidates(candidates);
00239 
00240   }
00241 }
00242 
00243 
00244 
00245 void CkfTrajectoryBuilder::updateTrajectory( TempTrajectory& traj,
00246                                              const TM& tm) const
00247 {
00248   TSOS predictedState = tm.predictedState();
00249   TM::ConstRecHitPointer hit = tm.recHit();
00250  
00251   if ( hit->isValid()) {
00252     TM tmp = TM( predictedState, theUpdator->update( predictedState, *hit),
00253                  hit, tm.estimate(), tm.layer()); 
00254     traj.push(tmp );
00255   }
00256   else {
00257     traj.push( TM( predictedState, hit, 0, tm.layer()));
00258   }
00259 }
00260 
00261 
00262 void 
00263 CkfTrajectoryBuilder::findCompatibleMeasurements( const TempTrajectory& traj, 
00264                                                   std::vector<TrajectoryMeasurement> & result) const
00265 {
00266   int invalidHits = 0;
00267   std::pair<TSOS,std::vector<const DetLayer*> > stateAndLayers = findStateAndLayers(traj);
00268   if (stateAndLayers.second.empty()) return;
00269 
00270   vector<const DetLayer*>::iterator layerBegin = stateAndLayers.second.begin();
00271   vector<const DetLayer*>::iterator layerEnd  = stateAndLayers.second.end();
00272   LogDebug("CkfPattern")<<"looping on "<< stateAndLayers.second.size()<<" layers.";
00273   for (vector<const DetLayer*>::iterator il = layerBegin; 
00274        il != layerEnd; il++) {
00275 
00276     LogDebug("CkfPattern")<<"looping on a layer in findCompatibleMeasurements.\n last layer: "<<traj.lastLayer()<<" current layer: "<<(*il);
00277 
00278     TSOS stateToUse = stateAndLayers.first;
00279     if ((*il)==traj.lastLayer())
00280       {
00281         LogDebug("CkfPattern")<<" self propagating in findCompatibleMeasurements.\n from: \n"<<stateToUse;
00282         //self navigation case
00283         // go to a middle point first
00284         TransverseImpactPointExtrapolator middle;
00285         GlobalPoint center(0,0,0);
00286         stateToUse = middle.extrapolate(stateToUse, center, *theForwardPropagator);
00287         
00288         if (!stateToUse.isValid()) continue;
00289         LogDebug("CkfPattern")<<"to: "<<stateToUse;
00290       }
00291     
00292     vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements((**il),stateToUse, *theForwardPropagator, *theEstimator);
00293     
00294     if ( !tmp.empty()) {
00295       if ( result.empty()) result = tmp;
00296       else {
00297         // keep one dummy TM at the end, skip the others
00298         result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
00299       }
00300       invalidHits++;
00301     }
00302   }
00303 
00304   // sort the final result, keep dummy measurements at the end
00305   if ( result.size() > 1) {
00306     sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
00307   }
00308 
00309   LogDebug("CkfPattern")<<"starting from:\n"
00310                         <<"x: "<<stateAndLayers.first.globalPosition()<<"\n"
00311                         <<"p: "<<stateAndLayers.first.globalMomentum()<<"\n"
00312                         <<PrintoutHelper::dumpMeasurements(result);
00313 
00314 #ifdef DEBUG_INVALID
00315   bool afterInvalid = false;
00316   for (vector<TM>::const_iterator i=result.begin();
00317        i!=result.end(); i++) {
00318     if ( ! i->recHit().isValid()) afterInvalid = true;
00319     if (afterInvalid && i->recHit().isValid()) {
00320       edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: valid hit after invalid!" ;
00321     }
00322   }
00323 #endif
00324 
00325   //analyseMeasurements( result, traj);
00326 
00327 }
00328