CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/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 void CkfTrajectoryBuilder::setEvent(const edm::Event& event) const
00060 {
00061   theMeasurementTracker->update(event);
00062 }
00063 
00064 CkfTrajectoryBuilder::TrajectoryContainer 
00065 CkfTrajectoryBuilder::trajectories(const TrajectorySeed& seed) const
00066 {  
00067   TrajectoryContainer result;
00068   result.reserve(5);
00069   trajectories(seed, result);
00070   return result;
00071 }
00072 
00073 /*
00074   void CkfTrajectoryBuilder::rememberSeedAndTrajectories(const TrajectorySeed& seed,
00075   CkfTrajectoryBuilder::TrajectoryContainer &result) const
00076   {
00077   
00078   //result ----> theCachedTrajectories
00079   //every first iteration on event. forget about everything that happened before
00080   if (edm::Service<UpdaterService>()->checkOnce(theUniqueName)) 
00081   theCachedTrajectories.clear();
00082   
00083   if (seed.nHits()==0) return;
00084   
00085   //then remember those trajectories
00086   for (TrajectoryContainer::iterator traj=result.begin();
00087   traj!=result.end(); ++traj) {
00088   theCachedTrajectories.insert(std::make_pair(seed.recHits().first->geographicalId(),*traj));
00089   }  
00090   }
00091   
00092   bool CkfTrajectoryBuilder::sharedSeed(const TrajectorySeed& s1,const TrajectorySeed& s2) const{
00093   //quit right away on nH=0
00094   if (s1.nHits()==0 || s2.nHits()==0) return false;
00095   //quit right away if not the same number of hits
00096   if (s1.nHits()!=s2.nHits()) return false;
00097   TrajectorySeed::range r1=s1.recHits();
00098   TrajectorySeed::range r2=s2.recHits();
00099   TrajectorySeed::const_iterator i1,i2;
00100   TrajectorySeed::const_iterator & i1_e=r1.second,&i2_e=r2.second;
00101   TrajectorySeed::const_iterator & i1_b=r1.first,&i2_b=r2.first;
00102   //quit right away if first detId does not match. front exist because of ==0 ->quit test
00103   if(i1_b->geographicalId() != i2_b->geographicalId()) return false;
00104   //then check hit by hit if they are the same
00105   for (i1=i1_b,i2=i2_b;i1!=i1_e,i2!=i2_e;++i1,++i2){
00106   if (!i1->sharesInput(&(*i2),TrackingRecHit::all)) return false;
00107   }
00108   return true;
00109   }
00110   bool CkfTrajectoryBuilder::seedAlreadyUsed(const TrajectorySeed& seed,
00111   CkfTrajectoryBuilder::TempTrajectoryContainer &candidates) const
00112   {
00113   //theCachedTrajectories ---> candidates
00114   if (seed.nHits()==0) return false;
00115   bool answer=false;
00116   pair<SharedTrajectory::const_iterator, SharedTrajectory::const_iterator> range = 
00117   theCachedTrajectories.equal_range(seed.recHits().first->geographicalId());
00118   SharedTrajectory::const_iterator trajP;
00119   for (trajP = range.first; trajP!=range.second;++trajP){
00120   //check whether seeds are identical     
00121   if (sharedSeed(trajP->second.seed(),seed)){
00122   candidates.push_back(trajP->second);
00123   answer=true;
00124   }//already existing trajectory shares the seed.   
00125   }//loop already made trajectories      
00126   
00127   return answer;
00128   }
00129 */
00130 
00131 void
00132 CkfTrajectoryBuilder::trajectories(const TrajectorySeed& seed, CkfTrajectoryBuilder::TrajectoryContainer &result) const
00133 {  
00134   // analyseSeed( seed);
00135   /*
00136     if (theSharedSeedCheck){
00137     TempTrajectoryContainer candidates;
00138     if (seedAlreadyUsed(seed,candidates))
00139     {
00140     //start with those candidates already made before
00141     limitedCandidates(candidates,result);
00142     //and quit
00143     return;
00144     }
00145     }
00146   */
00147 
00148   TempTrajectory startingTraj = createStartingTrajectory( seed );
00149 
00152   limitedCandidates( startingTraj, result);
00153 
00154   /*
00155     //and remember what you just did
00156     if (theSharedSeedCheck)  rememberSeedAndTrajectories(seed,result);
00157   */
00158 
00159   // analyseResult(result);
00160 }
00161 
00162 void CkfTrajectoryBuilder::
00163 limitedCandidates( TempTrajectory& startingTraj,
00164                    TrajectoryContainer& result) const
00165 {
00166   TempTrajectoryContainer candidates;
00167   candidates.push_back( startingTraj);
00168   limitedCandidates(candidates,result);
00169 }
00170 
00171 void CkfTrajectoryBuilder::
00172 limitedCandidates( TempTrajectoryContainer &candidates,
00173                    TrajectoryContainer& result) const
00174 {
00175   unsigned int nIter=1;
00176   //  TempTrajectoryContainer candidates; // = TrajectoryContainer();
00177   TempTrajectoryContainer newCand; // = TrajectoryContainer();
00178   //  candidates.push_back( startingTraj);
00179 
00180   while ( !candidates.empty()) {
00181 
00182     newCand.clear();
00183     for (TempTrajectoryContainer::iterator traj=candidates.begin();
00184          traj!=candidates.end(); traj++) {
00185       std::vector<TM> meas;
00186       findCompatibleMeasurements(*traj, meas);
00187 
00188       // --- method for debugging
00189       if(!analyzeMeasurementsDebugger(*traj,meas,
00190                                       theMeasurementTracker,
00191                                       theForwardPropagator,theEstimator,
00192                                       theTTRHBuilder)) return;
00193       // ---
00194 
00195       if ( meas.empty()) {
00196         if ( qualityFilter( *traj)) addToResult( *traj, result);
00197       }
00198       else {
00199         std::vector<TM>::const_iterator last;
00200         if ( theAlwaysUseInvalidHits) last = meas.end();
00201         else {
00202           if (meas.front().recHit()->isValid()) {
00203             last = find_if( meas.begin(), meas.end(), RecHitIsInvalid());
00204           }
00205           else last = meas.end();
00206         }
00207 
00208         for( std::vector<TM>::const_iterator itm = meas.begin(); 
00209              itm != last; itm++) {
00210           TempTrajectory newTraj = *traj;
00211           updateTrajectory( newTraj, *itm);
00212 
00213           if ( toBeContinued(newTraj)) {
00214             newCand.push_back(newTraj);
00215           }
00216           else {
00217             if ( qualityFilter(newTraj)) addToResult( newTraj, result);
00219           }
00220         }
00221       }
00222 
00223       if ((int)newCand.size() > theMaxCand) {
00224         sort( newCand.begin(), newCand.end(), TrajCandLess<TempTrajectory>(theLostHitPenalty));
00225         newCand.erase( newCand.begin()+theMaxCand, newCand.end());
00226       }
00227     }
00228 
00229     if (theIntermediateCleaning) IntermediateTrajectoryCleaner::clean(newCand);
00230 
00231     candidates.swap(newCand);
00232     
00233     LogDebug("CkfPattern") <<result.size()<<" candidates after "<<nIter++<<" CKF iteration: \n"
00234                            <<PrintoutHelper::dumpCandidates(result)
00235                            <<"\n "<<candidates.size()<<" running candidates are: \n"
00236                            <<PrintoutHelper::dumpCandidates(candidates);
00237 
00238   }
00239 }
00240 
00241 
00242 
00243 void CkfTrajectoryBuilder::updateTrajectory( TempTrajectory& traj,
00244                                              const TM& tm) const
00245 {
00246   TSOS predictedState = tm.predictedState();
00247   TM::ConstRecHitPointer hit = tm.recHit();
00248  
00249   if ( hit->isValid()) {
00250     TM tmp = TM( predictedState, theUpdator->update( predictedState, *hit),
00251                  hit, tm.estimate(), tm.layer()); 
00252     traj.push(tmp );
00253   }
00254   else {
00255     traj.push( TM( predictedState, hit, 0, tm.layer()));
00256   }
00257 }
00258 
00259 
00260 void 
00261 CkfTrajectoryBuilder::findCompatibleMeasurements( const TempTrajectory& traj, 
00262                                                   std::vector<TrajectoryMeasurement> & result) const
00263 {
00264   int invalidHits = 0;
00265   std::pair<TSOS,std::vector<const DetLayer*> > stateAndLayers = findStateAndLayers(traj);
00266   if (stateAndLayers.second.empty()) return;
00267 
00268   vector<const DetLayer*>::iterator layerBegin = stateAndLayers.second.begin();
00269   vector<const DetLayer*>::iterator layerEnd  = stateAndLayers.second.end();
00270   LogDebug("CkfPattern")<<"looping on "<< stateAndLayers.second.size()<<" layers.";
00271   for (vector<const DetLayer*>::iterator il = layerBegin; 
00272        il != layerEnd; il++) {
00273 
00274     LogDebug("CkfPattern")<<"looping on a layer in findCompatibleMeasurements.\n last layer: "<<traj.lastLayer()<<" current layer: "<<(*il);
00275 
00276     TSOS stateToUse = stateAndLayers.first;
00277     if ((*il)==traj.lastLayer())
00278       {
00279         LogDebug("CkfPattern")<<" self propagating in findCompatibleMeasurements.\n from: \n"<<stateToUse;
00280         //self navigation case
00281         // go to a middle point first
00282         TransverseImpactPointExtrapolator middle;
00283         GlobalPoint center(0,0,0);
00284         stateToUse = middle.extrapolate(stateToUse, center, *theForwardPropagator);
00285         
00286         if (!stateToUse.isValid()) continue;
00287         LogDebug("CkfPattern")<<"to: "<<stateToUse;
00288       }
00289     
00290     vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements((**il),stateToUse, *theForwardPropagator, *theEstimator);
00291     
00292     if ( !tmp.empty()) {
00293       if ( result.empty()) result = tmp;
00294       else {
00295         // keep one dummy TM at the end, skip the others
00296         result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
00297       }
00298       invalidHits++;
00299     }
00300   }
00301 
00302   // sort the final result, keep dummy measurements at the end
00303   if ( result.size() > 1) {
00304     sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
00305   }
00306 
00307   LogDebug("CkfPattern")<<"starting from:\n"
00308                         <<"x: "<<stateAndLayers.first.globalPosition()<<"\n"
00309                         <<"p: "<<stateAndLayers.first.globalMomentum()<<"\n"
00310                         <<PrintoutHelper::dumpMeasurements(result);
00311 
00312 #ifdef DEBUG_INVALID
00313   bool afterInvalid = false;
00314   for (vector<TM>::const_iterator i=result.begin();
00315        i!=result.end(); i++) {
00316     if ( ! i->recHit().isValid()) afterInvalid = true;
00317     if (afterInvalid && i->recHit().isValid()) {
00318       edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: valid hit after invalid!" ;
00319     }
00320   }
00321 #endif
00322 
00323   //analyseMeasurements( result, traj);
00324 
00325 }
00326