CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTracker/CkfPattern/src/GroupedCkfTrajectoryBuilder.cc

Go to the documentation of this file.
00001 
00002 #include "RecoTracker/CkfPattern/interface/GroupedCkfTrajectoryBuilder.h"
00003 #include "TrajectorySegmentBuilder.h"
00004 
00005 
00006 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00007 #include "TrackingTools/PatternTools/interface/TrajMeasLessEstim.h"
00008 
00009 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00010 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
00011 #include "TrackingTools/TrackFitters/interface/KFTrajectoryFitter.h"
00012 #include "RecoTracker/CkfPattern/interface/GroupedTrajCandLess.h"
00013 #include "TrackingTools/TrajectoryFiltering/interface/RegionalTrajectoryFilter.h"
00014 #include "TrackingTools/PatternTools/interface/TempTrajectory.h"
00015 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00016 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00017 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00018 #include "TrackingTools/DetLayers/interface/DetGroup.h"
00019 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00020 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00021 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00022 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00023 
00024 #include "TrackingTools/PatternTools/interface/TrajectoryStateUpdator.h"
00025 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00026 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00027 #include "TrackingTools/PatternTools/interface/TrajMeasLessEstim.h"
00028 #include "TrackingTools/TrajectoryState/interface/BasicSingleTrajectoryState.h"
00029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00030 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00031 
00032 // only included for RecHit comparison operator:
00033 #include "TrackingTools/TrajectoryCleaning/interface/TrajectoryCleanerBySharedHits.h"
00034 
00035 // for looper reconstruction
00036 #include "TrackingTools/GeomPropagators/interface/HelixBarrelCylinderCrossing.h"
00037 #include "TrackingTools/GeomPropagators/interface/HelixBarrelPlaneCrossingByCircle.h"
00038 #include "TrackingTools/GeomPropagators/interface/HelixArbitraryPlaneCrossing.h"
00039 
00040 
00041 #include <algorithm> 
00042 
00043 using namespace std;
00044 
00045 //#define DBG2_GCTB
00046 
00047 //#define STANDARD_INTERMEDIARYCLEAN
00048 
00049 #ifdef STANDARD_INTERMEDIARYCLEAN
00050 #include "RecoTracker/CkfPattern/interface/IntermediateTrajectoryCleaner.h"
00051 #endif
00052 
00053 /* ====== B.M. to be ported layer ===========
00054 #ifdef DBG_GCTB
00055 #include "RecoTracker/CkfPattern/src/ShowCand.h"
00056 #endif
00057 // #define DBG2_GCTB
00058 #ifdef DBG2_GCTB
00059 #include "RecoTracker/CkfPattern/src/SimIdPrinter.h"
00060 #include "Tracker/TkDebugTools/interface/LayerFinderByDet.h"
00061 #include "Tracker/TkLayout/interface/TkLayerName.h"
00062 #endif
00063 =================================== */
00064 
00065 
00066 GroupedCkfTrajectoryBuilder::
00067 GroupedCkfTrajectoryBuilder(const edm::ParameterSet&              conf,
00068                             const TrajectoryStateUpdator*         updator,
00069                             const Propagator*                     propagatorAlong,
00070                             const Propagator*                     propagatorOpposite,
00071                             const Chi2MeasurementEstimatorBase*   estimator,
00072                             const TransientTrackingRecHitBuilder* recHitBuilder,
00073                             const MeasurementTracker*             measurementTracker,
00074                             const TrajectoryFilter*               filter,
00075                             const TrajectoryFilter*               inOutFilter):
00076 
00077 
00078   BaseCkfTrajectoryBuilder(conf,
00079                            updator, propagatorAlong,propagatorOpposite,
00080                            estimator, recHitBuilder, measurementTracker, filter, inOutFilter)
00081 {
00082   // fill data members from parameters (eventually data members could be dropped)
00083   //
00084   theMaxCand                  = conf.getParameter<int>("maxCand");
00085 
00086   theLostHitPenalty           = conf.getParameter<double>("lostHitPenalty");
00087   theFoundHitBonus            = conf.getParameter<double>("foundHitBonus");
00088   theIntermediateCleaning     = conf.getParameter<bool>("intermediateCleaning");
00089   theAlwaysUseInvalid         = conf.getParameter<bool>("alwaysUseInvalidHits");
00090   theLockHits                 = conf.getParameter<bool>("lockHits");
00091   theBestHitOnly              = conf.getParameter<bool>("bestHitOnly");
00092   theMinNrOf2dHitsForRebuild  = 2;
00093   theRequireSeedHitsInRebuild = conf.getParameter<bool>("requireSeedHitsInRebuild");
00094   theMinNrOfHitsForRebuild    = max(0,conf.getParameter<int>("minNrOfHitsForRebuild"));
00095   maxPtForLooperReconstruction     = conf.existsAs<double>("maxPtForLooperReconstruction") ? 
00096     conf.getParameter<double>("maxPtForLooperReconstruction") : 0;
00097   maxDPhiForLooperReconstruction     = conf.existsAs<double>("maxDPhiForLooperReconstruction") ? 
00098     conf.getParameter<double>("maxDPhiForLooperReconstruction") : 2.0;
00099 
00100 
00101   /* ======= B.M. to be ported layer ===========
00102   bool setOK = thePropagator->setMaxDirectionChange(1.6);
00103   if (!setOK) 
00104     cout  << "GroupedCkfTrajectoryBuilder WARNING: "
00105           << "propagator does not support setMaxDirectionChange" 
00106           << endl;
00107   //   addStopCondition(theMinPtStopCondition);
00108 
00109   theConfigurableCondition = createAlgo<TrajectoryFilter>(componentConfig("StopCondition"));
00110   ===================================== */
00111 
00112 }
00113 
00114 /*
00115   void GroupedCkfTrajectoryBuilder::setEvent(const edm::Event& event) const
00116   {
00117   theMeasurementTracker->update(event);
00118 }
00119 */
00120 
00121 GroupedCkfTrajectoryBuilder::TrajectoryContainer 
00122 GroupedCkfTrajectoryBuilder::trajectories (const TrajectorySeed& seed) const 
00123 {
00124   TrajectoryContainer ret; 
00125   ret.reserve(10);
00126   buildTrajectories(seed, ret, 0);
00127   return ret; 
00128 }
00129 
00130 GroupedCkfTrajectoryBuilder::TrajectoryContainer 
00131 GroupedCkfTrajectoryBuilder::trajectories (const TrajectorySeed& seed, 
00132                                            const TrackingRegion& region) const
00133 {
00134   TrajectoryContainer ret; 
00135   ret.reserve(10);
00136   RegionalTrajectoryFilter regionalCondition(region);
00137   buildTrajectories(seed, ret, &regionalCondition);
00138   return ret; 
00139 }
00140 
00141 void 
00142 GroupedCkfTrajectoryBuilder::trajectories (const TrajectorySeed& seed, GroupedCkfTrajectoryBuilder::TrajectoryContainer &ret) const 
00143 {
00144   buildTrajectories(seed,ret,0);
00145 }
00146 
00147 void
00148 GroupedCkfTrajectoryBuilder::trajectories (const TrajectorySeed& seed, 
00149                                             GroupedCkfTrajectoryBuilder::TrajectoryContainer &ret,
00150                                             const TrackingRegion& region) const
00151 {
00152   RegionalTrajectoryFilter regionalCondition(region);
00153   buildTrajectories(seed,ret,&regionalCondition);
00154 }
00155 
00156 void  
00157 GroupedCkfTrajectoryBuilder::rebuildSeedingRegion(const TrajectorySeed& seed,
00158                                                   TrajectoryContainer& result) const {    
00159   TempTrajectory startingTraj = createStartingTrajectory( seed);
00160   TempTrajectoryContainer work;
00161 
00162   TrajectoryContainer final;
00163 
00164   work.reserve(result.size());
00165   for (TrajectoryContainer::iterator traj=result.begin();
00166        traj!=result.end(); ++traj) {
00167     if(traj->isValid()) work.push_back(TempTrajectory(*traj));
00168   }
00169 
00170   rebuildSeedingRegion(startingTraj,work);
00171   final.reserve(work.size());
00172 
00173   for (TempTrajectoryContainer::iterator traj=work.begin();
00174        traj!=work.end(); ++traj) {
00175     final.push_back(traj->toTrajectory());
00176   }
00177   
00178   result.swap(final);
00179   
00180 }
00181 
00182 void
00183 GroupedCkfTrajectoryBuilder::buildTrajectories (const TrajectorySeed& seed,
00184                                                 GroupedCkfTrajectoryBuilder::TrajectoryContainer &result,
00185                                                 const TrajectoryFilter* regionalCondition) const
00186 {
00187   //
00188   // Build trajectory outwards from seed
00189   //
00190 
00191   analyseSeed( seed);
00192 
00193   TempTrajectory startingTraj = createStartingTrajectory( seed);
00194 
00195   work_.clear();
00196   const bool inOut = true;
00197   groupedLimitedCandidates( startingTraj, regionalCondition, theForwardPropagator, inOut, work_);
00198   if ( work_.empty() )  return ;
00199 
00200 
00201 
00202   /*  rebuilding is de-coupled from standard building
00203   //
00204   // try to additional hits in the seeding region
00205   //
00206   if ( theMinNrOfHitsForRebuild>0 ) {
00207     // reverse direction
00208     //thePropagator->setPropagationDirection(oppositeDirection(seed.direction()));
00209     // rebuild part of the trajectory
00210     rebuildSeedingRegion(startingTraj,work);
00211   }
00212 
00213   */
00214 
00215   result.reserve(work_.size());
00216   for (TempTrajectoryContainer::const_iterator it = work_.begin(), ed = work_.end(); it != ed; ++it) {
00217       result.push_back( it->toTrajectory() );
00218   }
00219 
00220   work_.clear(); 
00221   if (work_.capacity() > work_MaxSize_) {  TempTrajectoryContainer().swap(work_); work_.reserve(work_MaxSize_/2); }
00222 
00223   analyseResult(result);
00224 
00225   LogDebug("CkfPattern")<< "GroupedCkfTrajectoryBuilder: returning result of size " << result.size();
00226 
00227 }
00228 
00229 
00230 void 
00231 GroupedCkfTrajectoryBuilder::groupedLimitedCandidates (TempTrajectory& startingTraj, 
00232                                                        const TrajectoryFilter* regionalCondition,
00233                                                        const Propagator* propagator, 
00234                                                        bool inOut,
00235                                                        TempTrajectoryContainer& result) const
00236 {
00237   unsigned int nIter=1;
00238   TempTrajectoryContainer candidates;
00239   TempTrajectoryContainer newCand;
00240   candidates.push_back( startingTraj);
00241 
00242   while ( !candidates.empty()) {
00243 
00244     newCand.clear();
00245     for (TempTrajectoryContainer::iterator traj=candidates.begin();
00246          traj!=candidates.end(); traj++) {
00247       if ( !advanceOneLayer(*traj, regionalCondition, propagator, inOut, newCand, result) ) {
00248         LogDebug("CkfPattern")<< "GCTB: terminating after advanceOneLayer==false";
00249         continue;
00250       }
00251 
00252       LogDebug("CkfPattern")<<"newCand(1): after advanced one layer:\n"<<PrintoutHelper::dumpCandidates(newCand);
00253 
00254       if ((int)newCand.size() > theMaxCand) {
00255         //ShowCand()(newCand);
00256 
00257         sort( newCand.begin(), newCand.end(), GroupedTrajCandLess(theLostHitPenalty,theFoundHitBonus));
00258         newCand.erase( newCand.begin()+theMaxCand, newCand.end());
00259       }
00260       LogDebug("CkfPattern")<<"newCand(2): after removing extra candidates.\n"<<PrintoutHelper::dumpCandidates(newCand);
00261     }
00262 
00263     LogDebug("CkfPattern") << "newCand.size() at end = " << newCand.size();
00264 /*
00265     if (theIntermediateCleaning) {
00266       candidates.clear();
00267       candidates = groupedIntermediaryClean(newCand);
00268     } else {
00269       candidates.swap(newCand);
00270     }
00271 */
00272     if (theIntermediateCleaning) {
00273 #ifdef STANDARD_INTERMEDIARYCLEAN
00274         IntermediateTrajectoryCleaner::clean(newCand);  
00275 #else 
00276         groupedIntermediaryClean(newCand);
00277 #endif  
00278 
00279     }   
00280     candidates.swap(newCand);
00281 
00282     LogDebug("CkfPattern") <<"candidates(3): "<<result.size()<<" candidates after "<<nIter++<<" groupedCKF iteration: \n"
00283                            <<PrintoutHelper::dumpCandidates(result)
00284                            <<"\n "<<candidates.size()<<" running candidates are: \n"
00285                            <<PrintoutHelper::dumpCandidates(candidates);
00286   }
00287 }
00288 
00289 std::string whatIsTheNextStep(TempTrajectory& traj , std::pair<TrajectoryStateOnSurface,std::vector<const DetLayer*> >& stateAndLayers){
00290   std::stringstream buffer;
00291   vector<const DetLayer*> & nl = stateAndLayers.second;
00292 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
00293 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
00294   //B.M. TkLayerName layerName;
00295   //B.M. buffer << "Started from " << layerName(traj.lastLayer()) 
00296   const BarrelDetLayer* sbdl = dynamic_cast<const BarrelDetLayer*>(traj.lastLayer());
00297   const ForwardDetLayer* sfdl = dynamic_cast<const ForwardDetLayer*>(traj.lastLayer());
00298   if (sbdl) {
00299     buffer << "Started from " << traj.lastLayer() << " r=" << sbdl->specificSurface().radius() 
00300            << " phi=" << sbdl->specificSurface().phi() << endl;
00301   } else if (sfdl) {
00302     buffer << "Started from " << traj.lastLayer() << " z " << sfdl->specificSurface().position().z()
00303            << " phi " << sfdl->specificSurface().phi() << endl;
00304   }
00305   buffer << "Trying to go to";
00306   for ( vector<const DetLayer*>::iterator il=nl.begin();
00307         il!=nl.end(); il++){ 
00308     //B.M. buffer << " " << layerName(*il)  << " " << *il << endl;
00309     const BarrelDetLayer* bdl = dynamic_cast<const BarrelDetLayer*>(*il);
00310     const ForwardDetLayer* fdl = dynamic_cast<const ForwardDetLayer*>(*il);
00311     
00312     if (bdl) buffer << " r " << bdl->specificSurface().radius() << endl;
00313     if (fdl) buffer << " z " << fdl->specificSurface().position().z() << endl;
00314     //buffer << " " << *il << endl;   
00315   }
00316   return buffer.str();
00317 }
00318 
00319 std::string whatIsTheStateToUse(TrajectoryStateOnSurface &initial, TrajectoryStateOnSurface & stateToUse, const DetLayer * l){
00320   std::stringstream buffer;
00321   buffer << "GCTB: starting from " 
00322          << " r / phi / z = " << stateToUse.globalPosition().perp()
00323          << " / " << stateToUse.globalPosition().phi()
00324          << " / " << stateToUse.globalPosition().z() 
00325          << " , pt / phi / pz /charge = " 
00326          << stateToUse.globalMomentum().perp() << " / "  
00327          << stateToUse.globalMomentum().phi() << " / " 
00328          << stateToUse.globalMomentum().z() << " / " 
00329          << stateToUse.charge()
00330          << " for layer at "<< l << endl;
00331   buffer << "     errors:";
00332   for ( int i=0; i<5; i++ )  buffer << " " << sqrt(stateToUse.curvilinearError().matrix()(i,i));
00333   buffer << endl;
00334   
00335   //buffer << "GCTB: starting from r / phi / z = " << initial.globalPosition().perp()
00336   //<< " / " << initial.globalPosition().phi()
00337   //<< " / " << initial.globalPosition().z() << " , pt / pz = " 
00338   //<< initial.globalMomentum().perp() << " / " 
00339   //<< initial.globalMomentum().z() << " for layer at "
00340   //<< l << endl;
00341   //buffer << "     errors:";
00342   //for ( int i=0; i<5; i++ )  buffer << " " << sqrt(initial.curvilinearError().matrix()(i,i));
00343   //buffer << endl;
00344   return buffer.str();
00345 }
00346 
00347 
00348 bool 
00349 GroupedCkfTrajectoryBuilder::advanceOneLayer (TempTrajectory& traj, 
00350                                               const TrajectoryFilter* regionalCondition, 
00351                                               const Propagator* propagator,
00352                                               bool inOut,
00353                                               TempTrajectoryContainer& newCand, 
00354                                               TempTrajectoryContainer& result) const
00355 {
00356   std::pair<TSOS,std::vector<const DetLayer*> > stateAndLayers = findStateAndLayers(traj);
00357   if(maxPtForLooperReconstruction>0){
00358     if(
00359        //stateAndLayers.second.size()==0 &&
00360        traj.lastLayer()->location()==0 && 
00361        stateAndLayers.first.globalMomentum().perp()<maxPtForLooperReconstruction &&
00362        stateAndLayers.first.globalMomentum().perp()>0.3)
00363       stateAndLayers.second.push_back(traj.lastLayer());
00364   }
00365 
00366   vector<const DetLayer*>::iterator layerBegin = stateAndLayers.second.begin();
00367   vector<const DetLayer*>::iterator layerEnd   = stateAndLayers.second.end();
00368 
00369   //   if (nl.empty()) {
00370   //     addToResult(traj,result,inOut);
00371   //     return false;
00372   //   }
00373   
00374   LogDebug("CkfPattern")<<whatIsTheNextStep(traj, stateAndLayers);
00375   
00376   bool foundSegments(false);
00377   bool foundNewCandidates(false);
00378   for ( vector<const DetLayer*>::iterator il=layerBegin; 
00379         il!=layerEnd; il++) {
00380 
00381     TSOS stateToUse = stateAndLayers.first;
00382 
00383     double dPhiCacheForLoopersReconstruction(0);
00384     if ((*il)==traj.lastLayer()){
00385       if(maxPtForLooperReconstruction>0){
00386         // ------ For loopers reconstruction
00387         //cout<<" self propagating in advanceOneLayer (for loopers) \n";
00388         const BarrelDetLayer* sbdl = dynamic_cast<const BarrelDetLayer*>(traj.lastLayer());
00389         if(sbdl){
00390           HelixBarrelCylinderCrossing cylinderCrossing(stateToUse.globalPosition(),
00391                                                        stateToUse.globalMomentum(),
00392                                                        stateToUse.transverseCurvature(),
00393                                                        propagator->propagationDirection(),
00394                                                        sbdl->specificSurface());
00395           if(!cylinderCrossing.hasSolution()) continue;
00396           GlobalPoint starting = stateToUse.globalPosition();
00397           GlobalPoint target1 = cylinderCrossing.position1();
00398           GlobalPoint target2 = cylinderCrossing.position2();
00399           
00400           GlobalPoint farther = fabs(starting.phi()-target1.phi()) > fabs(starting.phi()-target2.phi()) ?
00401             target1 : target2;
00402         
00403           const Bounds& bounds( sbdl->specificSurface().bounds());
00404           float length = bounds.length() / 2.f;
00405         
00406           /*
00407             cout << "starting: " << starting << endl;
00408             cout << "target1: " << target1 << endl;
00409             cout << "target2: " << target2 << endl;
00410             cout << "dphi: " << (target1.phi()-target2.phi()) << endl;
00411             cout << "length: " << length << endl;
00412           */
00413         
00414           /*
00415           float deltaZ = bounds.thickness()/2.f/fabs(tan(stateToUse.globalDirection().theta()) ) ;
00416           if(stateToUse.hasError())
00417             deltaZ += 3*sqrt(stateToUse.cartesianError().position().czz());
00418           if( fabs(farther.z()) > length + deltaZ ) continue;
00419           */
00420           if(fabs(farther.z())*0.95>length) continue;
00421 
00422           Geom::Phi<double> tmpDphi = target1.phi()-target2.phi();
00423           if(fabs(tmpDphi)>maxDPhiForLooperReconstruction) continue;
00424           GlobalPoint target((target1.x()+target2.x())/2,
00425                              (target1.y()+target2.y())/2,
00426                              (target1.z()+target2.z())/2);
00427           
00428           //cout << "target: " << target << endl;
00429           
00430 
00431           
00432           TransverseImpactPointExtrapolator extrapolator;
00433           stateToUse = extrapolator.extrapolate(stateToUse, target, *propagator);
00434           if (!stateToUse.isValid()) continue; //SK: consider trying the original? probably not
00435 
00436           //dPhiCacheForLoopersReconstruction = fabs(target1.phi()-target2.phi())*2.;
00437           dPhiCacheForLoopersReconstruction = fabs(tmpDphi);
00438           traj.incrementLoops();
00439         }else{
00440           continue;
00441         }
00442       }else{
00443         // ------ For cosmics reconstruction
00444         LogDebug("CkfPattern")<<" self propagating in advanceOneLayer.\n from: \n"<<stateToUse;
00445         //self navigation case
00446         // go to a middle point first
00447         TransverseImpactPointExtrapolator middle;
00448         GlobalPoint center(0,0,0);
00449         stateToUse = middle.extrapolate(stateToUse, center, *theForwardPropagator);
00450         
00451         if (!stateToUse.isValid()) continue;
00452         LogDebug("CkfPattern")<<"to: "<<stateToUse;
00453       }
00454     }
00455     
00456     unsigned int maxCandidates = theMaxCand > 21 ? theMaxCand*2 : 42; //limit the number of returned segments
00457     TrajectorySegmentBuilder layerBuilder(theMeasurementTracker,
00458                                           theLayerMeasurements,
00459                                           **il,*propagator,
00460                                           *theUpdator,*theEstimator,
00461                                           theLockHits,theBestHitOnly, maxCandidates);
00462 
00463     LogDebug("CkfPattern")<<whatIsTheStateToUse(stateAndLayers.first,stateToUse,*il);
00464     
00465     TempTrajectoryContainer segments=
00466       layerBuilder.segments(stateToUse);
00467 
00468     LogDebug("CkfPattern")<< "GCTB: number of segments = " << segments.size();
00469 
00470     if ( !segments.empty() )  foundSegments = true;
00471 
00472     for ( TempTrajectoryContainer::const_iterator is=segments.begin();
00473           is!=segments.end(); is++ ) {
00474       //
00475       // assume "invalid hit only" segment is last in list
00476       //
00477       const TempTrajectory::DataContainer & measurements = is->measurements();
00478       if ( !theAlwaysUseInvalid && is!=segments.begin() && measurements.size()==1 && 
00479            (measurements.front().recHit()->getType() == TrackingRecHit::missing) )  break;
00480       //
00481       // create new candidate
00482       //
00483       TempTrajectory newTraj(traj);
00484       traj.setDPhiCacheForLoopersReconstruction(dPhiCacheForLoopersReconstruction);
00485       
00486       //----  avoid to add the same hits more than once in the trajectory ----
00487       bool toBeRejected(false);
00488       for(const TempTrajectory::DataContainer::const_iterator revIt = measurements.rbegin(); 
00489           revIt!=measurements.rend(); --revIt){
00490         int tmpCounter(0);
00491         for(const TempTrajectory::DataContainer::const_iterator newTrajMeasIt = newTraj.measurements().rbegin(); 
00492             newTrajMeasIt != newTraj.measurements().rend(); --newTrajMeasIt){
00493           //if(tmpCounter==2) break;
00494           if(revIt->recHit()->geographicalId()==newTrajMeasIt->recHit()->geographicalId()){
00495             toBeRejected=true;
00496             break;
00497           }
00498           tmpCounter++;
00499         }
00500       }
00501 
00502       if(!toBeRejected){
00503         //newTraj.push(*is);
00504         //std::cout << "DEBUG: newTraj after push found,lost: " 
00505         //        << newTraj.foundHits() << " , " 
00506         //        << newTraj.lostHits() << " , "
00507         //        << newTraj.measurements().size() << std::endl;
00508       }else{    
00509         /*cout << "WARNING: neglect candidate because it contains the same hit twice \n";
00510         cout << "-- discarded track's pt,eta,#found: " 
00511              << newTraj.lastMeasurement().updatedState().globalMomentum().perp() << " , "
00512              << newTraj.lastMeasurement().updatedState().globalMomentum().eta() << " , "
00513              << newTraj.foundHits() << "\n";
00514         */
00515         continue; //Are we sure about this????
00516       }
00517       // ------------------------
00518 
00519 
00520       newTraj.push(*is);
00521       //GIO// for ( vector<TM>::const_iterator im=measurements.begin();
00522       //GIO//        im!=measurements.end(); im++ )  newTraj.push(*im);
00523       //if ( toBeContinued(newTraj,regionalCondition) ) { TOBE FIXED
00524       if ( toBeContinued(newTraj, inOut) ) {
00525         // Have added one more hit to track candidate
00526         
00527         LogDebug("CkfPattern")<<"GCTB: adding updated trajectory to candidates: inOut="<<inOut<<" hits="<<newTraj.foundHits();
00528 
00529         newCand.push_back(newTraj);
00530         foundNewCandidates = true;
00531       }
00532       else {
00533         // Have finished building this track. Check if it passes cuts.
00534 
00535         LogDebug("CkfPattern")<< "GCTB: adding completed trajectory to results if passes cuts: inOut="<<inOut<<" hits="<<newTraj.foundHits();
00536 
00537         addToResult(newTraj, result, inOut);
00538       }
00539     }
00540   }
00541 
00542   if ( !foundSegments ){
00543     LogDebug("CkfPattern")<< "GCTB: adding input trajectory to result";
00544     addToResult(traj, result, inOut);
00545   }
00546   return foundNewCandidates;
00547 }
00548 
00549 //TempTrajectoryContainer
00550 void
00551 GroupedCkfTrajectoryBuilder::groupedIntermediaryClean (TempTrajectoryContainer& theTrajectories) const 
00552 {
00553   //if (theTrajectories.empty()) return TrajectoryContainer();
00554   //TrajectoryContainer result;
00555   if (theTrajectories.empty()) return;  
00556   //RecHitEqualByChannels recHitEqualByChannels(false, false);
00557   int firstLayerSize, secondLayerSize;
00558   vector<const DetLayer*> firstLayers, secondLayers;
00559 
00560   for (TempTrajectoryContainer::iterator firstTraj=theTrajectories.begin();
00561        firstTraj!=(theTrajectories.end()-1); firstTraj++) {
00562 
00563     if ( (!firstTraj->isValid()) ||
00564          (!firstTraj->lastMeasurement().recHit()->isValid()) ) continue;
00565     const TempTrajectory::DataContainer & firstMeasurements = firstTraj->measurements();
00566     layers(firstMeasurements, firstLayers);
00567     firstLayerSize = firstLayers.size();
00568     if ( firstLayerSize<4 )  continue;
00569 
00570     for (TempTrajectoryContainer::iterator secondTraj=(firstTraj+1);
00571        secondTraj!=theTrajectories.end(); secondTraj++) {
00572 
00573       if ( (!secondTraj->isValid()) ||
00574            (!secondTraj->lastMeasurement().recHit()->isValid()) ) continue;
00575       const TempTrajectory::DataContainer & secondMeasurements = secondTraj->measurements();
00576       layers(secondMeasurements, secondLayers);
00577       secondLayerSize = secondLayers.size();
00578       //
00579       // only candidates using the same last 3 layers are compared
00580       //
00581       if ( firstLayerSize!=secondLayerSize )  continue;
00582       if ( firstLayers[0]!=secondLayers[0] ||
00583            firstLayers[1]!=secondLayers[1] ||
00584            firstLayers[2]!=secondLayers[2] )  continue;
00585 
00586       TempTrajectory::DataContainer::const_iterator im1 = firstMeasurements.rbegin();
00587       TempTrajectory::DataContainer::const_iterator im2 = secondMeasurements.rbegin();
00588       //
00589       // check for identical hits in the last layer
00590       //
00591       bool unequal(false);
00592       const DetLayer* layerPtr = firstLayers[0];
00593       while ( im1!=firstMeasurements.rend()&&im2!=secondMeasurements.rend() ) {
00594         if ( im1->layer()!=layerPtr || im2->layer()!=layerPtr )  break;
00595         if ( !(im1->recHit()->isValid()) || !(im2->recHit()->isValid()) ||
00596              !im1->recHit()->hit()->sharesInput(im2->recHit()->hit(), TrackingRecHit::some) ) {
00598           unequal = true;
00599           break;
00600         }
00601         --im1;
00602         --im2;
00603       }
00604       if ( im1==firstMeasurements.rend() || im2==secondMeasurements.rend() ||
00605            im1->layer()==layerPtr || im2->layer()==layerPtr || unequal )  continue;
00606       //
00607       // check for invalid hits in the layer -2
00608       // compare only candidates with invalid / valid combination
00609       //
00610       layerPtr = firstLayers[1];
00611       bool firstValid(true);
00612       while ( im1!=firstMeasurements.rend()&&im1->layer()==layerPtr ) {
00613         if ( !im1->recHit()->isValid() )  firstValid = false;
00614         --im1;
00615       }
00616       bool secondValid(true);
00617       while ( im2!=secondMeasurements.rend()&&im2->layer()==layerPtr ) {
00618         if ( !im2->recHit()->isValid() )  secondValid = false;
00619         --im2;
00620       }
00621       if ( !tkxor(firstValid,secondValid) )  continue;
00622       //
00623       // ask for identical hits in layer -3
00624       //
00625       unequal = false;
00626       layerPtr = firstLayers[2];
00627       while ( im1!=firstMeasurements.rend()&&im2!=secondMeasurements.rend() ) {
00628         if ( im1->layer()!=layerPtr || im2->layer()!=layerPtr )  break;
00629         if ( !(im1->recHit()->isValid()) || !(im2->recHit()->isValid()) ||
00630              !im1->recHit()->hit()->sharesInput(im2->recHit()->hit(), TrackingRecHit::some) ) {
00632           unequal = true;
00633           break;
00634         }
00635         --im1;
00636         --im2;
00637       }
00638       if ( im1==firstMeasurements.rend() || im2==secondMeasurements.rend() ||
00639            im1->layer()==layerPtr || im2->layer()==layerPtr || unequal )  continue;
00640 
00641       if ( !firstValid ) {
00642         firstTraj->invalidate();
00643         break;
00644       }
00645       else {
00646         secondTraj->invalidate();
00647         break;
00648       }
00649     }
00650   }
00651 /*
00652   for (TempTrajectoryContainer::const_iterator it = theTrajectories.begin();
00653        it != theTrajectories.end(); it++) {
00654     if(it->isValid()) result.push_back( *it);
00655   }
00656 
00657   return result;
00658 */
00659   theTrajectories.erase(std::remove_if( theTrajectories.begin(),theTrajectories.end(),
00660                                         std::not1(std::mem_fun_ref(&TempTrajectory::isValid))),
00661  //                                     boost::bind(&TempTrajectory::isValid,_1)), 
00662                         theTrajectories.end());
00663 }
00664 
00665 void
00666 GroupedCkfTrajectoryBuilder::layers (const TempTrajectory::DataContainer& measurements,
00667                                      vector<const DetLayer*> &result) const 
00668 {
00669   result.clear();
00670 
00671   if ( measurements.empty() )  return ;
00672 
00673   result.push_back(measurements.back().layer());
00674   TempTrajectory::DataContainer::const_iterator ifirst = measurements.rbegin();
00675   --ifirst;      
00676   for ( TempTrajectory::DataContainer::const_iterator im=ifirst;
00677         im!=measurements.rend(); --im ) {
00678     if ( im->layer()!=result.back() )  result.push_back(im->layer());
00679   }
00680 
00681   for (vector<const DetLayer*>::const_iterator iter = result.begin(); iter != result.end(); iter++){
00682     if (!*iter) edm::LogWarning("CkfPattern")<< "Warning: null det layer!! ";
00683   }
00684 }
00685 
00686 void
00687 GroupedCkfTrajectoryBuilder::rebuildSeedingRegion(TempTrajectory& startingTraj,
00688                                                   TempTrajectoryContainer& result) const
00689 {
00690   //
00691   // Rebuilding of trajectories. Candidates are taken from result,
00692   // which will be replaced with the solutions after rebuild
00693   // (assume vector::swap is more efficient than building new container)
00694   //
00695   LogDebug("CkfPattern")<< "Starting to rebuild " << result.size() << " tracks";
00696   //
00697   // Fitter (need to create it here since the propagation direction
00698   // might change between different starting trajectories)
00699   //
00700   KFTrajectoryFitter fitter(&(*theBackwardPropagator),&updator(),&estimator());
00701   //
00702   TempTrajectoryContainer reFitted;
00703   TrajectorySeed::range rseedHits = startingTraj.seed().recHits();
00704   std::vector<const TrackingRecHit*> seedHits;
00705   //seedHits.insert(seedHits.end(), rseedHits.first, rseedHits.second);
00706   //for (TrajectorySeed::recHitContainer::const_iterator iter = rseedHits.first; iter != rseedHits.second; iter++){
00707   //    seedHits.push_back(&*iter);
00708   //}
00709 
00710   //unsigned int nSeed(seedHits.size());
00711   unsigned int nSeed(rseedHits.second-rseedHits.first);
00712   //seedHits.reserve(nSeed);
00713   TempTrajectoryContainer rebuiltTrajectories;
00714   for ( TempTrajectoryContainer::iterator it=result.begin();
00715         it!=result.end(); it++ ) {
00716     //
00717     // skip candidates which are not exceeding the seed size 
00718     // (e.g. because no Tracker layers outside seeding region) 
00719     //
00720 
00721     if ( it->measurements().size()<=startingTraj.measurements().size() ) {
00722       rebuiltTrajectories.push_back(*it);
00723       LogDebug("CkfPattern")<< "RebuildSeedingRegion skipped as in-out trajectory does not exceed seed size.";
00724       continue;
00725     }
00726     //
00727     // Refit - keep existing trajectory in case fit is not possible
00728     // or fails
00729     //
00730 
00731     backwardFit(*it,nSeed,fitter,reFitted,seedHits);
00732     if ( reFitted.size()!=1 ) {
00733       rebuiltTrajectories.push_back(*it);
00734       LogDebug("CkfPattern")<< "RebuildSeedingRegion skipped as backward fit failed";
00735       //                            << "after reFitted.size() " << reFitted.size();
00736       continue;
00737     }
00738     //LogDebug("CkfPattern")<<"after reFitted.size() " << reFitted.size();
00739     //
00740     // Rebuild seeding part. In case it fails: keep initial trajectory
00741     // (better to drop it??)
00742     //
00743     int nRebuilt =
00744       rebuildSeedingRegion (seedHits,reFitted.front(),rebuiltTrajectories);
00745 
00746     if ( nRebuilt==0 ) it->invalidate();  // won't use original in-out track
00747 
00748     if ( nRebuilt<=0 ) rebuiltTrajectories.push_back(*it);
00749 
00750   }
00751   //
00752   // Replace input trajectories with new ones
00753   //
00754   result.swap(rebuiltTrajectories);
00755   result.erase(std::remove_if( result.begin(),result.end(),
00756                                std::not1(std::mem_fun_ref(&TempTrajectory::isValid))),
00757                result.end());
00758 }
00759 
00760 int
00761 GroupedCkfTrajectoryBuilder::rebuildSeedingRegion(const std::vector<const TrackingRecHit*>& seedHits, 
00762                                                   TempTrajectory& candidate,
00763                                                   TempTrajectoryContainer& result) const 
00764 {
00765   //
00766   // Starting from track found by in-out tracking phase, extrapolate it inwards through
00767   // the seeding region if possible in towards smaller Tracker radii, searching for additional
00768   // hits.
00769   // The resulting trajectories are returned in result,
00770   // the count is the return value.
00771   //
00772   TempTrajectoryContainer rebuiltTrajectories;
00773 #ifdef DBG2_GCTB
00774 /*  const LayerFinderByDet layerFinder;
00775   if ( !seedHits.empty() && seedHits.front().isValid() ) {
00776     DetLayer* seedLayer = layerFinder(seedHits.front().det());
00777     cout << "Seed hit at " << seedHits.front().globalPosition()
00778          << " " << seedLayer << endl;
00779     cout << "Started from " 
00780          << candidate.lastMeasurement().updatedState().globalPosition().perp() << " "
00781          << candidate.lastMeasurement().updatedState().globalPosition().z() << endl;
00782     pair<bool,TrajectoryStateOnSurface> layerComp(false,TrajectoryStateOnSurface());
00783     if ( seedLayer ) layerComp =
00784       seedLayer->compatible(candidate.lastMeasurement().updatedState(),
00785                               propagator(),estimator());
00786     pair<bool,TrajectoryStateOnSurface> detComp =
00787       seedHits.front().det().compatible(candidate.lastMeasurement().updatedState(),
00788                                         propagator(),estimator());
00789     cout << "  layer compatibility = " << layerComp.first;
00790     cout << "  det compatibility = " << detComp.first;
00791     if ( detComp.first ) {
00792       cout << "  estimate = " 
00793            << estimator().estimate(detComp.second,seedHits.front()).second ;
00794     }
00795     cout << endl;
00796   }*/
00797   cout << "Before backward building: #measurements = " 
00798        << candidate.measurements().size() ; //<< endl;;
00799 #endif
00800   //
00801   // Use standard building with standard cuts. Maybe better to use different
00802   // cuts from "forward" building (e.g. no check on nr. of invalid hits)?
00803   //
00804   const bool inOut = false;
00805   groupedLimitedCandidates(candidate, (const TrajectoryFilter*)0, theBackwardPropagator, inOut, rebuiltTrajectories);
00806 
00807   LogDebug("CkfPattern")<<" After backward building: "<<PrintoutHelper::dumpCandidates(rebuiltTrajectories);
00808 
00809   //
00810   // Check & count resulting candidates
00811   //
00812   int nrOfTrajectories(0);
00813   bool orig_ok = false;
00814   //const RecHitEqualByChannels recHitEqual(false,false);
00815   //vector<TM> oldMeasurements(candidate.measurements());
00816   for ( TempTrajectoryContainer::iterator it=rebuiltTrajectories.begin();
00817         it!=rebuiltTrajectories.end(); it++ ) {
00818 
00819     TempTrajectory::DataContainer newMeasurements(it->measurements());
00820     //
00821     // Verify presence of seeding hits?
00822     //
00823     if ( theRequireSeedHitsInRebuild ) {
00824       orig_ok = true;
00825       // no hits found (and possibly some invalid hits discarded): drop track
00826       if ( newMeasurements.size()<=candidate.measurements().size() ){  
00827         LogDebug("CkfPattern") << "newMeasurements.size()<=candidate.measurements().size()";
00828         continue;
00829       } 
00830       // verify presence of hits
00831       //GIO//if ( !verifyHits(newMeasurements.begin()+candidate.measurements().size(),
00832       //GIO//                  newMeasurements.end(),seedHits) ){
00833       if ( !verifyHits(newMeasurements.rbegin(), 
00834                        newMeasurements.size() - candidate.measurements().size(),
00835                        seedHits) ){
00836         LogDebug("CkfPattern")<< "seed hits not found in rebuild";
00837         continue; 
00838       }
00839     }
00840     //
00841     // construct final trajectory in the right order
00842     //
00843     TempTrajectory reversedTrajectory(it->seed(),it->seed().direction());
00844     reversedTrajectory.setNLoops(it->nLoops());
00845     for (TempTrajectory::DataContainer::const_iterator im=newMeasurements.rbegin(), ed = newMeasurements.rend();
00846           im != ed; --im ) {
00847       reversedTrajectory.push(*im);
00848     }
00849     // save & count result
00850     result.push_back(reversedTrajectory);
00851     nrOfTrajectories++;
00852 
00853     LogDebug("CkgPattern")<<"New traj direction = " << reversedTrajectory.direction()<<"\n"
00854                           <<PrintoutHelper::dumpMeasurements(reversedTrajectory.measurements());
00855   }
00856   // If nrOfTrajectories = 0 and orig_ok = true, this means that a track was actually found on the
00857   // out-in step (meeting those requirements) but did not have the seed hits in it.
00858   // In this case when we return we will go ahead and use the original in-out track.
00859 
00860   // If nrOfTrajectories = 0 and orig_ok = false, this means that the out-in step failed to
00861   // find any track.  Two cases are a technical failure in fitting the original seed hits or
00862   // because the track did not meet the out-in criteria (which may be stronger than the out-in
00863   // criteria).  In this case we will NOT allow the original in-out track to be used.
00864 
00865   if ( (nrOfTrajectories == 0) && orig_ok ) {
00866     nrOfTrajectories = -1;
00867   }
00868   return nrOfTrajectories;
00869 }
00870 
00871 void
00872 GroupedCkfTrajectoryBuilder::backwardFit (TempTrajectory& candidate, unsigned int nSeed,
00873                                           const TrajectoryFitter& fitter,
00874                                           TempTrajectoryContainer& fittedTracks,
00875                                           std::vector<const TrackingRecHit*>& remainingHits) const
00876 {
00877   //
00878   // clear array of non-fitted hits
00879   //
00880   remainingHits.clear();
00881   fittedTracks.clear();
00882   //
00883   // skip candidates which are not exceeding the seed size
00884   // (e.g. Because no Tracker layers exist outside seeding region)
00885   //
00886   if ( candidate.measurements().size()<=nSeed ) {
00887     fittedTracks.clear();
00888     return;
00889   }
00890 
00891   LogDebug("CkfPattern")<<"nSeed " << nSeed << endl
00892                         << "Old traj direction = " << candidate.direction() << endl
00893                         <<PrintoutHelper::dumpMeasurements(candidate.measurements());
00894 
00895   //
00896   // backward fit trajectory.
00897   // (Will try to fit only hits outside the seeding region. However,
00898   // if there are not enough of these, it will also use the seeding hits).
00899   //
00900   TempTrajectory::DataContainer oldMeasurements(candidate.measurements());
00901 //   int nOld(oldMeasurements.size());
00902 //   const unsigned int nHitAllMin(5);
00903 //   const unsigned int nHit2dMin(2);
00904   unsigned int nHit(0);    // number of valid hits after seeding region
00905   //unsigned int nHit2d(0);  // number of valid hits after seeding region with 2D info
00906   // use all hits except the first n (from seed), but require minimum
00907   // specified in configuration.
00908   //  Swapped over next two lines.
00909   unsigned int nHitMin = max(candidate.foundHits()-nSeed,theMinNrOfHitsForRebuild);
00910   //  unsigned int nHitMin = oldMeasurements.size()-nSeed;
00911   // we want to rebuild only if the number of VALID measurements excluding the seed measurements is higher than the cut
00912   if (nHitMin<theMinNrOfHitsForRebuild){
00913         fittedTracks.clear();
00914         return;
00915   }
00916 
00917   LogDebug("CkfPattern")/* << "nHitMin " << nHitMin*/ <<"Sizes: " << oldMeasurements.size() << " / ";
00918   //
00919   // create input trajectory for backward fit
00920   //
00921   Trajectory fwdTraj(candidate.seed(),oppositeDirection(candidate.direction()));
00922   fwdTraj.setNLoops(candidate.nLoops());
00923   //const TrajectorySeed seed = TrajectorySeed(PTrajectoryStateOnDet(), TrajectorySeed::recHitContainer(), oppositeDirection(candidate.direction()));
00924   //Trajectory fwdTraj(seed, oppositeDirection(candidate.direction()));
00925   std::vector<const DetLayer*> bwdDetLayer; 
00926   for ( TempTrajectory::DataContainer::const_iterator im=oldMeasurements.rbegin();
00927         im!=oldMeasurements.rend(); --im) {
00928     const TrackingRecHit* hit = im->recHit()->hit();
00929     //
00930     // add hits until required number is reached
00931     //
00932     if ( nHit<nHitMin ){//|| nHit2d<theMinNrOf2dHitsForRebuild ) {
00933       fwdTraj.push(*im);
00934       bwdDetLayer.push_back(im->layer());
00935       //
00936       // count valid / 2D hits
00937       //
00938       if ( hit->isValid() ) {
00939         nHit++;
00940         //if ( hit.isMatched() ||
00941         //     hit.det().detUnits().front()->type().module()==pixel )
00942         //nHit2d++;
00943       }
00944     }
00945     //if (nHit==nHitMin) lastBwdDetLayer=im->layer();   
00946     //
00947     // keep remaining (valid) hits for verification
00948     //
00949     else if ( hit->isValid() ) {
00950       //std::cout << "Adding a remaining hit" << std::endl;
00951       remainingHits.push_back(hit);
00952     }
00953   }
00954   //
00955   // Fit only if required number of valid hits can be used
00956   //
00957   if ( nHit<nHitMin ){  //|| nHit2d<theMinNrOf2dHitsForRebuild ) {
00958     fittedTracks.clear();
00959     return;
00960   }
00961   //
00962   // Do the backward fit (important: start from scaled, not random cov. matrix!)
00963   //
00964   TrajectoryStateOnSurface firstTsos(fwdTraj.firstMeasurement().updatedState());
00965   //cout << "firstTsos "<< firstTsos << endl;
00966   firstTsos.rescaleError(10.);
00967   //TrajectoryContainer bwdFitted(fitter.fit(fwdTraj.seed(),fwdTraj.recHits(),firstTsos));
00968   TrajectoryContainer bwdFitted(fitter.fit(
00969                 TrajectorySeed(PTrajectoryStateOnDet(), TrajectorySeed::recHitContainer(), oppositeDirection(candidate.direction())),
00970                 fwdTraj.recHits(),firstTsos));
00971   if (bwdFitted.size()){
00972     LogDebug("CkfPattern")<<"Obtained " << bwdFitted.size() << " bwdFitted trajectories with measurement size " << bwdFitted.front().measurements().size();
00973         TempTrajectory fitted(fwdTraj.seed(), fwdTraj.direction());
00974         fitted.setNLoops(fwdTraj.nLoops());
00975         vector<TM> tmsbf = bwdFitted.front().measurements();
00976         int iDetLayer=0;
00977         //this is ugly but the TM in the fitted track do not contain the DetLayer.
00978         //So we have to cache the detLayer pointers and replug them in.
00979         //For the backward building it would be enaugh to cache the last DetLayer, 
00980         //but for the intermediary cleaning we need all
00981         for ( vector<TM>::const_iterator im=tmsbf.begin();im!=tmsbf.end(); im++ ) {
00982                 fitted.push(TM( (*im).forwardPredictedState(),
00983                                 (*im).backwardPredictedState(),
00984                                 (*im).updatedState(),
00985                                 (*im).recHit(),
00986                                 (*im).estimate(),
00987                                 bwdDetLayer[iDetLayer]));
00988 
00989                 LogDebug("CkfPattern")<<PrintoutHelper::dumpMeasurement(*im);
00990                 iDetLayer++;
00991         }
00992 /*
00993         TM lastMeas = bwdFitted.front().lastMeasurement();
00994         fitted.pop();
00995         fitted.push(TM(lastMeas.forwardPredictedState(), 
00996                                lastMeas.backwardPredictedState(), 
00997                                lastMeas.updatedState(),
00998                                lastMeas.recHit(),
00999                                lastMeas.estimate(),
01000                                lastBwdDetLayer));*/
01001         fittedTracks.push_back(fitted);
01002   }
01003   //
01004   // save result
01005   //
01006   //fittedTracks.swap(bwdFitted);
01007   //cout << "Obtained " << fittedTracks.size() << " fittedTracks trajectories with measurement size " << fittedTracks.front().measurements().size() << endl;
01008 }
01009 
01010 bool
01011 GroupedCkfTrajectoryBuilder::verifyHits (TempTrajectory::DataContainer::const_iterator rbegin,
01012                                          size_t maxDepth,
01013                                          const std::vector<const TrackingRecHit*>& hits) const
01014 {
01015   //
01016   // verify presence of the seeding hits
01017   //
01018   LogDebug("CkfPattern")<<"Checking for " << hits.size() << " hits in "
01019                         << maxDepth << " measurements" << endl;
01020 
01021   TempTrajectory::DataContainer::const_iterator rend = rbegin; 
01022   while (maxDepth > 0) { --maxDepth; --rend; }
01023   for ( vector<const TrackingRecHit*>::const_iterator ir=hits.begin();
01024         ir!=hits.end(); ir++ ) {
01025     // assume that all seeding hits are valid!
01026     bool foundHit(false);
01027     for ( TempTrajectory::DataContainer::const_iterator im=rbegin; im!=rend; --im ) {
01028       if ( im->recHit()->isValid() && (*ir)->sharesInput(im->recHit()->hit(), TrackingRecHit::some) ) {
01029         foundHit = true;
01030         break;
01031       }
01032     }
01033     if ( !foundHit )  return false;
01034   }
01035   return true;
01036 }
01037 
01038 
01039 
01040