CMS 3D CMS Logo

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