CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/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   buildTrajectories(seed, result,nullptr);
00151 }
00152 
00153 TempTrajectory CkfTrajectoryBuilder::buildTrajectories (const TrajectorySeed&seed,
00154                                                         TrajectoryContainer &result,
00155                                                         const TrajectoryFilter*) const {
00156   
00157   TempTrajectory startingTraj = createStartingTrajectory( seed );
00158   
00161   limitedCandidates(seed, startingTraj, result);
00162   
00163   return startingTraj;
00164 
00165   /*
00166   //and remember what you just did
00167   if (theSharedSeedCheck)  rememberSeedAndTrajectories(seed,result);
00168   */
00169   
00170   // analyseResult(result);
00171 }
00172 
00173 void CkfTrajectoryBuilder::
00174 limitedCandidates(const TrajectorySeed&seed, TempTrajectory& startingTraj,
00175                    TrajectoryContainer& result) const
00176 {
00177   TempTrajectoryContainer candidates;
00178   candidates.push_back( startingTraj);
00179   boost::shared_ptr<const TrajectorySeed>  sharedSeed(new TrajectorySeed(seed));
00180   limitedCandidates(sharedSeed, candidates,result);
00181 }
00182 
00183 void CkfTrajectoryBuilder::
00184 limitedCandidates(const boost::shared_ptr<const TrajectorySeed> & sharedSeed, TempTrajectoryContainer &candidates,
00185                    TrajectoryContainer& result) const
00186 {
00187   unsigned int nIter=1;
00188   TempTrajectoryContainer newCand; // = TrajectoryContainer();
00189 
00190  
00191   while ( !candidates.empty()) {
00192 
00193     newCand.clear();
00194     for (TempTrajectoryContainer::iterator traj=candidates.begin();
00195          traj!=candidates.end(); traj++) {
00196       std::vector<TM> meas;
00197       findCompatibleMeasurements(*sharedSeed, *traj, meas);
00198 
00199       // --- method for debugging
00200       if(!analyzeMeasurementsDebugger(*traj,meas,
00201                                       theMeasurementTracker,
00202                                       theForwardPropagator,theEstimator,
00203                                       theTTRHBuilder)) return;
00204       // ---
00205 
00206       if ( meas.empty()) {
00207         if ( qualityFilter( *traj)) addToResult(sharedSeed, *traj, result);
00208       }
00209       else {
00210         std::vector<TM>::const_iterator last;
00211         if ( theAlwaysUseInvalidHits) last = meas.end();
00212         else {
00213           if (meas.front().recHit()->isValid()) {
00214             last = find_if( meas.begin(), meas.end(), RecHitIsInvalid());
00215           }
00216           else last = meas.end();
00217         }
00218 
00219         for( std::vector<TM>::const_iterator itm = meas.begin(); 
00220              itm != last; itm++) {
00221           TempTrajectory newTraj = *traj;
00222           updateTrajectory( newTraj, *itm);
00223 
00224           if ( toBeContinued(newTraj)) {
00225             newCand.push_back(newTraj);
00226           }
00227           else {
00228             if ( qualityFilter(newTraj)) addToResult(sharedSeed, newTraj, result);
00230           }
00231         }
00232       }
00233 
00234       if ((int)newCand.size() > theMaxCand) {
00235         sort( newCand.begin(), newCand.end(), TrajCandLess<TempTrajectory>(theLostHitPenalty));
00236         newCand.erase( newCand.begin()+theMaxCand, newCand.end());
00237       }
00238     }
00239 
00240     if (theIntermediateCleaning) IntermediateTrajectoryCleaner::clean(newCand);
00241 
00242     candidates.swap(newCand);
00243     
00244     LogDebug("CkfPattern") <<result.size()<<" candidates after "<<nIter++<<" CKF iteration: \n"
00245                            <<PrintoutHelper::dumpCandidates(result)
00246                            <<"\n "<<candidates.size()<<" running candidates are: \n"
00247                            <<PrintoutHelper::dumpCandidates(candidates);
00248 
00249   }
00250 }
00251 
00252 
00253 
00254 void CkfTrajectoryBuilder::updateTrajectory( TempTrajectory& traj,
00255                                              const TM& tm) const
00256 {
00257   TSOS predictedState = tm.predictedState();
00258   TM::ConstRecHitPointer hit = tm.recHit();
00259  
00260   if ( hit->isValid()) {
00261     TM tmp = TM( predictedState, theUpdator->update( predictedState, *hit),
00262                  hit, tm.estimate(), tm.layer()); 
00263     traj.push(tmp );
00264   }
00265   else {
00266     traj.push( TM( predictedState, hit, 0, tm.layer()));
00267   }
00268 }
00269 
00270 
00271 void 
00272 CkfTrajectoryBuilder::findCompatibleMeasurements(const TrajectorySeed&seed,
00273                                                  const TempTrajectory& traj, 
00274                                                  std::vector<TrajectoryMeasurement> & result) const
00275 {
00276   int invalidHits = 0;
00277   std::pair<TSOS,std::vector<const DetLayer*> > stateAndLayers = findStateAndLayers(traj);
00278   if (stateAndLayers.second.empty()) return;
00279 
00280   vector<const DetLayer*>::iterator layerBegin = stateAndLayers.second.begin();
00281   vector<const DetLayer*>::iterator layerEnd  = stateAndLayers.second.end();
00282   LogDebug("CkfPattern")<<"looping on "<< stateAndLayers.second.size()<<" layers.";
00283   for (vector<const DetLayer*>::iterator il = layerBegin; 
00284        il != layerEnd; il++) {
00285 
00286     LogDebug("CkfPattern")<<"looping on a layer in findCompatibleMeasurements.\n last layer: "<<traj.lastLayer()<<" current layer: "<<(*il);
00287 
00288     TSOS stateToUse = stateAndLayers.first;
00289     if ((*il)==traj.lastLayer())
00290       {
00291         LogDebug("CkfPattern")<<" self propagating in findCompatibleMeasurements.\n from: \n"<<stateToUse;
00292         //self navigation case
00293         // go to a middle point first
00294         TransverseImpactPointExtrapolator middle;
00295         GlobalPoint center(0,0,0);
00296         stateToUse = middle.extrapolate(stateToUse, center, *theForwardPropagator);
00297         
00298         if (!stateToUse.isValid()) continue;
00299         LogDebug("CkfPattern")<<"to: "<<stateToUse;
00300       }
00301     
00302     vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements((**il),stateToUse, *theForwardPropagator, *theEstimator);
00303     
00304     if ( !tmp.empty()) {
00305       if ( result.empty()) result = tmp;
00306       else {
00307         // keep one dummy TM at the end, skip the others
00308         result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
00309       }
00310       invalidHits++;
00311     }
00312   }
00313 
00314   // sort the final result, keep dummy measurements at the end
00315   if ( result.size() > 1) {
00316     sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
00317   }
00318 
00319   LogDebug("CkfPattern")<<"starting from:\n"
00320                         <<"x: "<<stateAndLayers.first.globalPosition()<<"\n"
00321                         <<"p: "<<stateAndLayers.first.globalMomentum()<<"\n"
00322                         <<PrintoutHelper::dumpMeasurements(result);
00323 
00324 #ifdef DEBUG_INVALID
00325   bool afterInvalid = false;
00326   for (vector<TM>::const_iterator i=result.begin();
00327        i!=result.end(); i++) {
00328     if ( ! i->recHit().isValid()) afterInvalid = true;
00329     if (afterInvalid && i->recHit().isValid()) {
00330       edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: valid hit after invalid!" ;
00331     }
00332   }
00333 #endif
00334 
00335   //analyseMeasurements( result, traj);
00336 
00337 }
00338