CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/RecoTracker/SingleTrackPattern/src/CRackTrajectoryBuilder.cc

Go to the documentation of this file.
00001 //
00002 // Package:         RecoTracker/SingleTrackPattern
00003 // Class:           CRackTrajectoryBuilder
00004 // Original Author:  Michele Pioppi-INFN perugia
00005 //
00006 // Package:         RecoTracker/SingleTrackPattern
00007 // Class:           CRackTrajectoryBuilder
00008 // Original Author:  Michele Pioppi-INFN perugia
00009 
00010 #include <vector>
00011 #include <iostream>
00012 #include <cmath>
00013 
00014 #include "RecoTracker/SingleTrackPattern/interface/CRackTrajectoryBuilder.h"
00015 #include "DataFormats/Common/interface/OwnVector.h"
00016 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00017 #include "TrackingTools/TrajectoryState/interface/BasicSingleTrajectoryState.h"
00018 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h" 
00023 #include "TrackingTools/Records/interface/TransientRecHitRecord.h" 
00024 #include "TrackingTools/TrackFitters/interface/TrajectoryStateWithArbitraryError.h"
00025 
00026 
00027 //#include "RecoTracker/CkfPattern/interface/TransientInitialStateEstimator.h"
00028 
00029 using namespace std;
00030 
00031 
00032 CRackTrajectoryBuilder::CRackTrajectoryBuilder(const edm::ParameterSet& conf) : conf_(conf) { 
00033   //minimum number of hits per tracks
00034 
00035   theMinHits=conf_.getParameter<int>("MinHits");
00036   //cut on chi2
00037   chi2cut=conf_.getParameter<double>("Chi2Cut");
00038   edm::LogInfo("CosmicTrackFinder")<<"Minimum number of hits "<<theMinHits<<" Cut on Chi2= "<<chi2cut;
00039 
00040   debug_info=conf_.getUntrackedParameter<bool>("debug", false);
00041   fastPropagation=conf_.getUntrackedParameter<bool>("fastPropagation", false);
00042   useMatchedHits=conf_.getUntrackedParameter<bool>("useMatchedHits", true);
00043 
00044 
00045   geometry=conf_.getUntrackedParameter<std::string>("GeometricStructure","STANDARD");
00046 
00047   
00048 
00049 }
00050 
00051 
00052 CRackTrajectoryBuilder::~CRackTrajectoryBuilder() {
00053 //  delete theInitialState;
00054 }
00055 
00056 
00057 void CRackTrajectoryBuilder::init(const edm::EventSetup& es, bool seedplus){
00058 
00059 
00060 //  edm::ParameterSet tise_params = conf_.getParameter<edm::ParameterSet>("TransientInitialStateEstimatorParameters") ;
00061 // theInitialState          = new TransientInitialStateEstimator( es,tise_params);
00062 
00063   //services
00064   es.get<IdealMagneticFieldRecord>().get(magfield);
00065   es.get<TrackerDigiGeometryRecord>().get(tracker);
00066   
00067  
00068   
00069   if (seedplus) {        
00070     seed_plus=true;      
00071     thePropagator=      new PropagatorWithMaterial(alongMomentum,0.1057,&(*magfield) );          
00072     thePropagatorOp=    new PropagatorWithMaterial(oppositeToMomentum,0.1057,&(*magfield) );}    
00073   else {
00074     seed_plus=false;
00075     thePropagator=      new PropagatorWithMaterial(oppositeToMomentum,0.1057,&(*magfield) );    
00076     thePropagatorOp=    new PropagatorWithMaterial(alongMomentum,0.1057,&(*magfield) );
00077   }
00078   
00079   theUpdator=       new KFUpdator();
00080 //  theUpdator=       new KFStripUpdator();
00081   theEstimator=     new Chi2MeasurementEstimator(chi2cut);
00082   
00083 
00084   edm::ESHandle<TransientTrackingRecHitBuilder> theBuilder;
00085   std::string builderName = conf_.getParameter<std::string>("TTRHBuilder");   
00086   es.get<TransientRecHitRecord>().get(builderName,theBuilder);
00087   
00088 
00089   RHBuilder=   theBuilder.product();
00090 
00091 
00092 
00093 
00094   theFitter=        new KFTrajectoryFitter(*thePropagator,
00095                                            *theUpdator, 
00096                                            *theEstimator) ;
00097   
00098 
00099   theSmoother=      new KFTrajectorySmoother(*thePropagatorOp,
00100                                              *theUpdator,       
00101                                              *theEstimator);
00102   
00103 }
00104 
00105 
00106 void CRackTrajectoryBuilder::run(const TrajectorySeedCollection &collseed,
00107                                   const SiStripRecHit2DCollection &collstereo,
00108                                   const SiStripRecHit2DCollection &collrphi ,
00109                                   const SiStripMatchedRecHit2DCollection &collmatched,
00110                                   const SiPixelRecHitCollection &collpixel,
00111                                   const edm::EventSetup& es,
00112                                   edm::Event& e,
00113                                   vector<Trajectory> &trajoutput)
00114 {
00115 
00116   std::vector<Trajectory> trajSmooth;
00117   std::vector<Trajectory>::iterator trajIter;
00118   
00119   TrajectorySeedCollection::const_iterator iseed;
00120   unsigned int IS=0;
00121   for(iseed=collseed.begin();iseed!=collseed.end();iseed++){
00122     bool seedplus=((*iseed).direction()==alongMomentum);
00123     init(es,seedplus);
00124     hits.clear();
00125     trajFit.clear();
00126     trajSmooth.clear();
00127 
00128     Trajectory startingTraj = createStartingTrajectory(*iseed);
00129     Trajectory trajTmp; // used for the opposite direction
00130     Trajectory traj = startingTraj;
00131 
00132    
00133     //first propagate track in opposite direction ...
00134     seed_plus = !seed_plus;
00135     vector<const TrackingRecHit*> allHitsOppsite = SortHits(collstereo,collrphi,collmatched,collpixel, *iseed, true);
00136     seed_plus = !seed_plus;
00137     if (allHitsOppsite.size())
00138       {
00139         //there are hits which are above the seed,
00140         //cout << "Number of hits higher than seed " <<allHitsOppsite.size() << endl;   
00141         
00142         AddHit( traj, allHitsOppsite, thePropagatorOp);
00143         
00144         
00145         if (debug_info)
00146           {
00147             cout << "Hits in opposite direction..." << endl;
00148             TransientTrackingRecHit::RecHitContainer::const_iterator iHit;
00149             for ( iHit = hits.begin(); iHit!=hits.end(); iHit++)
00150               {
00151                 cout <<  (**iHit).globalPosition() << endl;
00152               }
00153           }
00154         
00155         
00156         // now add hist opposite to seed
00157         //now crate new starting trajectory
00158         reverse(hits.begin(), hits.end());
00159         if (thePropagatorOp->propagate(traj.firstMeasurement().updatedState(),
00160                                        tracker->idToDet((hits.front())->geographicalId())->surface()).isValid()){
00161           TSOS startingStateNew =  //TrajectoryStateWithArbitraryError()//
00162             (thePropagatorOp->propagate(traj.firstMeasurement().updatedState(),
00163                                         tracker->idToDet((hits.front())->geographicalId())->surface()));
00164           
00165           if (debug_info)
00166             {
00167               cout << "Hits in opposite direction reversed..." << endl;
00168               TransientTrackingRecHit::RecHitContainer::const_iterator iHit;
00169               for ( iHit = hits.begin(); iHit!=hits.end(); iHit++)
00170                 {
00171                   cout <<  (**iHit).globalPosition() << endl;
00172                 }
00173             }
00174           
00175           const TrajectorySeed& tmpseed=traj.seed();
00176           trajTmp = theFitter->fit(tmpseed,hits, startingStateNew ).front();
00177           
00178           if(debug_info){
00179             cout << "Debugging show fitted hits" << endl;           
00180                 std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> > hitsFit= trajTmp.recHits();
00181                 std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> >::const_iterator hit;
00182                 for(hit=hitsFit.begin();hit!=hitsFit.end();hit++){
00183                   
00184                   cout << RHBuilder->build( &(*(*hit)->hit()) )->globalPosition() << endl;
00185                 }
00186           }
00187           
00188         }
00189       }
00190     else 
00191       {
00192         if(debug_info) cout << "There are no hits in opposite direction ..." << endl;
00193       }
00194 
00195     vector<const TrackingRecHit*> allHits;
00196     if (trajTmp.foundHits())
00197       {
00198         traj = trajTmp;
00199         allHits= SortHits(collstereo,collrphi,collmatched,collpixel,*iseed, false);       
00200       }
00201     else
00202       {
00203         traj = startingTraj;
00204         hits.clear();
00205         allHits= SortHits(collstereo,collrphi,collmatched,collpixel,*iseed, true);        
00206       }
00207     
00208     AddHit( traj,allHits, thePropagator);
00209     
00210     
00211     if(debug_info){
00212       cout << "Debugging show All fitted hits" << endl;     
00213       std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> > hits= traj.recHits();
00214             std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> >::const_iterator hit;
00215             for(hit=hits.begin();hit!=hits.end();hit++){
00216               
00217               cout << RHBuilder->build( &(*(*hit)->hit()) )->globalPosition() << endl;
00218             }
00219 
00220             cout << qualityFilter( traj) << " <- quality filter good?" << endl;
00221           }
00222 
00223 
00224     if (debug_info) cout << "now do quality checks" << endl;
00225     if ( qualityFilter( traj) ) {
00226       const TrajectorySeed& tmpseed=traj.seed();
00227       std::pair<TrajectoryStateOnSurface, const GeomDet*> initState =  innerState(traj); //theInitialState->innerState(traj);
00228         if (initState.first.isValid())
00229               trajFit = theFitter->fit(tmpseed,hits, initState.first);
00230     }
00231 
00232     for (trajIter=trajFit.begin(); trajIter!=trajFit.end();trajIter++){
00233       trajSmooth=theSmoother->trajectories((*trajIter));
00234     }
00235     for (trajIter= trajSmooth.begin(); trajIter!=trajSmooth.end();trajIter++){
00236       if((*trajIter).isValid()){
00237 
00238         if (debug_info) cout << "adding track ... " << endl;
00239         trajoutput.push_back((*trajIter));
00240       }
00241     }
00242     delete theUpdator;
00243     delete theEstimator;
00244     delete thePropagator;
00245     delete thePropagatorOp;
00246     delete theFitter;
00247     delete theSmoother;
00248     //Only the first 30 seeds are considered
00249     if (IS>30) return;
00250     IS++;
00251 
00252   }
00253 }
00254 
00255 Trajectory CRackTrajectoryBuilder::createStartingTrajectory( const TrajectorySeed& seed) const
00256 {
00257   Trajectory result( seed, seed.direction());
00258   std::vector<TM> seedMeas = seedMeasurements(seed);
00259   if ( !seedMeas.empty()) {
00260     for (std::vector<TM>::const_iterator i=seedMeas.begin(); i!=seedMeas.end(); i++){
00261       result.push(*i);
00262     }
00263   }
00264  
00265   return result;
00266 }
00267 
00268 
00269 std::vector<TrajectoryMeasurement> 
00270 CRackTrajectoryBuilder::seedMeasurements(const TrajectorySeed& seed) const
00271 {
00272   std::vector<TrajectoryMeasurement> result;
00273   TrajectorySeed::range hitRange = seed.recHits();
00274   for (TrajectorySeed::const_iterator ihit = hitRange.first; 
00275        ihit != hitRange.second; ihit++) {
00276     //RC TransientTrackingRecHit* recHit = RHBuilder->build(&(*ihit));
00277     TransientTrackingRecHit::RecHitPointer recHit = RHBuilder->build(&(*ihit));
00278     const GeomDet* hitGeomDet = (&(*tracker))->idToDet( ihit->geographicalId());
00279     TSOS invalidState( new BasicSingleTrajectoryState( hitGeomDet->surface()));
00280 
00281     if (ihit == hitRange.second - 1) {
00282       TSOS  updatedState=startingTSOS(seed);
00283       result.push_back(TM( invalidState, updatedState, recHit));
00284 
00285     } 
00286     else {
00287       result.push_back(TM( invalidState, recHit));
00288     }
00289     
00290   }
00291 
00292   return result;
00293 }
00294 
00295 vector<const TrackingRecHit*> 
00296 CRackTrajectoryBuilder::SortHits(const SiStripRecHit2DCollection &collstereo,
00297                                   const SiStripRecHit2DCollection &collrphi ,
00298                                   const SiStripMatchedRecHit2DCollection &collmatched,
00299                                   const SiPixelRecHitCollection &collpixel,
00300                                   const TrajectorySeed &seed,
00301                                   const bool bAddSeedHits
00302                                   ){
00303 
00304 
00305   //The Hits with global y more than the seed are discarded
00306   //The Hits correspondign to the seed are discarded
00307   //At the end all the hits are sorted in y
00308   vector<const TrackingRecHit*> allHits;
00309 
00310   SiStripRecHit2DCollection::DataContainer::const_iterator istrip;
00311   TrajectorySeed::range hRange= seed.recHits();
00312   TrajectorySeed::const_iterator ihit;
00313   float yref=0.;
00314 
00315   if (debug_info) cout << "SEED " << startingTSOS(seed) << endl; 
00316   if (debug_info) cout << "seed hits size " << seed.nHits() << endl;
00317 
00318 
00319   //seed debugging:
00320  // GeomDet *seedDet =  tracker->idToDet(seed.startingState().detId());
00321 
00322 
00323 //  edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "SEED " << seed.startingState(); 
00324 
00325   edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "SEED " << startingTSOS(seed); 
00326 //  edm::LogInfo("CRackTrajectoryBuilder::SortHits" << "seed hits size " << seed.nHits();
00327 
00328 
00329 
00330   //  if (seed.nHits()<2)
00331   //  return allHits;
00332 
00333   float_t yMin=0.;
00334   float_t yMax=0.;
00335 
00336   int seedHitSize= hRange.second - hRange.first;
00337 
00338   vector <int> detIDSeedMatched (seedHitSize);
00339   vector <int> detIDSeedRphi (seedHitSize);
00340   vector <int> detIDSeedStereo (seedHitSize);
00341 
00342   for (ihit = hRange.first; 
00343        ihit != hRange.second; ihit++) {
00344 
00345      // need to find track with lowest (seed_plus)/ highest y (seed_minus)
00346     // split matched hits ...
00347    const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>(&(*ihit));
00348 
00349    yref=RHBuilder->build(&(*ihit))->globalPosition().y();
00350    if (ihit == hRange.first)
00351      {
00352         yMin = yref;
00353         yMax = yref; 
00354      }
00355 
00356    if (matchedhit)
00357      {
00358        auto m = matchedhit->monoHit();
00359        auto s = matchedhit->stereoHit();
00360        float_t yGlobRPhi   = RHBuilder->build(&m)->globalPosition().y();
00361        float_t yGlobStereo = RHBuilder->build(&s)->globalPosition().y();
00362 
00363        if (debug_info) cout << "Rphi ..." << yGlobRPhi << endl;
00364        if (debug_info) cout << "Stereo ..." << yGlobStereo << endl;
00365 
00366        if ( yGlobStereo < yMin ) yMin = yGlobStereo;
00367        if ( yGlobRPhi   < yMin ) yMin = yGlobRPhi;
00368 
00369        if ( yGlobStereo > yMax ) yMax = yGlobStereo;
00370        if ( yGlobRPhi   > yMax ) yMax = yGlobRPhi;
00371 
00372        detIDSeedMatched.push_back (  matchedhit->geographicalId().rawId() );
00373        detIDSeedRphi.push_back ( m.geographicalId().rawId() );
00374        detIDSeedStereo.push_back ( s.geographicalId().rawId() );
00375 
00376       if (bAddSeedHits)
00377         {
00378           if (useMatchedHits) 
00379             {
00380               
00381               hits.push_back((RHBuilder->build(&(*ihit)))); 
00382             }
00383           else
00384             {
00385               if ( ( (yGlobRPhi > yGlobStereo ) && seed_plus ) || ((yGlobRPhi < yGlobStereo ) && !seed_plus ))  
00386                 {
00387                   hits.push_back((RHBuilder->build(&m)));
00388                   hits.push_back((RHBuilder->build(&s)));
00389                 }
00390               else
00391                 {
00392                   hits.push_back((RHBuilder->build(&s)));     
00393                   hits.push_back((RHBuilder->build(&m)));
00394                   
00395                 }
00396             }
00397         }
00398      }
00399    else if (bAddSeedHits)
00400      {
00401        hits.push_back((RHBuilder->build(&(*ihit)))); 
00402        detIDSeedRphi.push_back ( ihit->geographicalId().rawId() );
00403        detIDSeedMatched.push_back ( -1 );
00404        detIDSeedStereo.push_back ( -1 );
00405        
00406      }
00407        
00408    if ( yref < yMin ) yMin = yref;
00409    if ( yref > yMax ) yMax = yref;
00410    
00411 //    if (bAddSeedHits)
00412 //      hits.push_back((RHBuilder->build(&(*ihit)))); 
00413 
00414     LogDebug("CosmicTrackFinder")<<"SEED HITS"<<RHBuilder->build(&(*ihit))->globalPosition();
00415     if (debug_info) cout <<"SEED HITS"<<RHBuilder->build(&(*ihit))->globalPosition() << endl;
00416 //    if (debug_info) cout <<"SEED HITS"<< seed.startingState().parameters() << endl;
00417 
00418 
00419   }
00420   
00421   yref = (seed_plus) ? yMin : yMax;
00422   
00423   if ((&collpixel)!=0){
00424     SiPixelRecHitCollection::DataContainer::const_iterator ipix;
00425     for(ipix=collpixel.data().begin();ipix!=collpixel.data().end();ipix++){
00426       float ych= RHBuilder->build(&(*ipix))->globalPosition().y();
00427       if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
00428         allHits.push_back(&(*ipix));
00429     }
00430   } 
00431   
00432 
00433   if (useMatchedHits) // use matched
00434     {
00435         //add the matched hits ...
00436       SiStripMatchedRecHit2DCollection::DataContainer::const_iterator istripm;
00437 
00438       if ((&collmatched)!=0){
00439         for(istripm=collmatched.data().begin();istripm!=collmatched.data().end();istripm++){
00440           float ych= RHBuilder->build(&(*istripm))->globalPosition().y();
00441 
00442           int cDetId=istripm->geographicalId().rawId();
00443           bool noSeedDet = ( detIDSeedMatched.end() == find (detIDSeedMatched.begin(), detIDSeedMatched.end(), cDetId ) ) ;
00444 
00445           if ( noSeedDet )
00446           if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
00447             {
00448           //if (debug_info) cout << "adding matched hit " << &(*istripm) << endl; 
00449               allHits.push_back(&(*istripm));
00450         }
00451         }
00452       }
00453 
00454    //add the rpi hits, but only accept hits that are not matched hits
00455   if ((&collrphi)!=0){
00456     for(istrip=collrphi.data().begin();istrip!=collrphi.data().end();istrip++){
00457       float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
00458       StripSubdetector monoDetId(istrip->geographicalId());
00459       if (monoDetId.partnerDetId())
00460         {
00461           edm::LogInfo("CRackTrajectoryBuilder::SortHits")  << "this det belongs to a glued det " << ych << endl;
00462           continue;
00463         }
00464           int cDetId=istrip->geographicalId().rawId();
00465           bool noSeedDet = ( detIDSeedRphi.end()== find (detIDSeedRphi.begin(), detIDSeedRphi.end(), cDetId ) ) ;
00466           if (noSeedDet)
00467       if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
00468         {
00469          
00470           bool hitIsUnique = true;
00471           //now 
00472            if ((&collmatched)!=0)
00473              for(istripm=collmatched.data().begin();istripm!=collmatched.data().end();istripm++)
00474                {
00475                  //              if ( isDifferentStripReHit2D ( *istrip, (istripm->stereoHit() ) ) == false)
00476                    if ( isDifferentStripReHit2D ( *istrip, (istripm->monoHit() ) ) == false)
00477                    {
00478                      hitIsUnique = false;
00479                      edm::LogInfo("CRackTrajectoryBuilder::SortHits")  << "rphi hit is in matched hits; y: " << ych << endl;
00480                      break;
00481                    }
00482                } //end loop over all matched
00483            if (hitIsUnique)
00484              {
00485                  //      if (debug_info) cout << "adding rphi hit " << &(*istrip) << endl; 
00486                 allHits.push_back(&(*istrip));   
00487              }
00488         }
00489     }
00490   }
00491 
00492   
00493   //add the stereo hits except the hits that are in the matched collection
00494   //update do not use unmatched rphi hist due to limitation of alignment framework
00495   //if (!useMatchedHits)
00496   //if ((&collstereo)!=0){
00497   //  for(istrip=collstereo.data().begin();istrip!=collstereo.data().end();istrip++){
00498   //    float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
00499   //
00500   //
00501   //      int cDetId = istrip->geographicalId().rawId();
00502   //        bool noSeedDet = ( detIDSeedStereo.end()== find (detIDSeedStereo.begin(), detIDSeedStereo.end(), cDetId ) ) ;
00503   //
00504   //        if (noSeedDet)
00505   //    if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
00506   //      {
00507   //
00508   //        bool hitIsUnique = true;
00509   //        //now 
00510   //         if ((&collmatched)!=0)
00511   //           for(istripm=collmatched.data().begin();istripm!=collmatched.data().end();istripm++)
00512   //             {
00513   //             if ( isDifferentStripReHit2D ( *istrip, * (istripm->stereoHit() ) ) == false)
00514   //               {
00515   //                 hitIsUnique = false;
00516   //                 edm::LogInfo("CRackTrajectoryBuilder::SortHits") << "stereo hit already in matched hits; y:  " << ych << endl;
00517   //                 break;
00518   //               }
00519   //             } //end loop over all stereo
00520   //         if (hitIsUnique)
00521   //           {
00522   //             
00523   //          //   if (debug_info) cout << "now I am adding a stero hit, either noise or not in overlap ...!!!!" << endl;
00524   //                allHits.push_back(&(*istrip));   
00525   //           }
00526   //      }
00527   //  }
00528   //}
00529     }
00530   else // dont use matched ...
00531     {
00532       
00533       if ((&collrphi)!=0){
00534         for(istrip=collrphi.data().begin();istrip!=collrphi.data().end();istrip++){
00535           float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
00536           if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
00537             allHits.push_back(&(*istrip));   
00538         }
00539       }
00540 
00541 
00542       if ((&collstereo)!=0){
00543         for(istrip=collstereo.data().begin();istrip!=collstereo.data().end();istrip++){
00544           float ych= RHBuilder->build(&(*istrip))->globalPosition().y();
00545           if ((seed_plus && (ych<yref)) || (!(seed_plus) && (ych>yref)))
00546             allHits.push_back(&(*istrip));
00547         }
00548       }
00549 
00550     }
00551 
00552 
00553 //  if (seed_plus){
00554 //    stable_sort(allHits.begin(),allHits.end(),CompareDetY_plus(*tracker));
00555 //  }
00556 //  else {
00557 //    stable_sort(allHits.begin(),allHits.end(),CompareDetY_minus(*tracker));
00558 //  }
00559 
00560 
00561   if (seed_plus){
00562       stable_sort(allHits.begin(),allHits.end(),CompareHitY_plus(*tracker));
00563   }
00564   else {
00565       stable_sort(allHits.begin(),allHits.end(),CompareHitY(*tracker));
00566   }
00567 
00568 
00569 
00570   if (debug_info) 
00571     {
00572       if (debug_info) cout << "all hits" << endl;
00573 
00574       //starting trajectory
00575       Trajectory startingTraj = createStartingTrajectory(seed);
00576 
00577       if (debug_info) cout << "START " << startingTraj.lastMeasurement().updatedState() << endl;
00578       if (debug_info) cout << "START  Err" << startingTraj.lastMeasurement().updatedState().localError().matrix() << endl;
00579 
00580 
00581       vector<const TrackingRecHit*>::iterator iHit;
00582       for (iHit=allHits.begin(); iHit<allHits.end(); iHit++)
00583         {
00584               GlobalPoint gphit=RHBuilder->build(*iHit)->globalPosition();
00585               if (debug_info) cout << "GH " << gphit << endl;
00586 
00587               //      tracker->idToDet((*iHit)->geographicalId())->surface();
00588 
00589       TSOS prSt = thePropagator->propagate(startingTraj.lastMeasurement().updatedState(),
00590                           tracker->idToDet((*iHit)->geographicalId())->surface());
00591               
00592               if(prSt.isValid()) 
00593                 {
00594       if (debug_info) cout << "PR " << prSt.globalPosition() << endl;
00595       //if (debug_info) cout << "PR  Err" << prSt.localError().matrix() << endl;
00596                   
00597                 }
00598               else
00599                 {
00600                   if (debug_info) cout << "not valid" << endl;
00601                 }
00602         }
00603               if (debug_info) cout << "all hits end" << endl;
00604     }
00605 
00606 
00607   return allHits;
00608 }
00609 
00610 TrajectoryStateOnSurface
00611 CRackTrajectoryBuilder::startingTSOS(const TrajectorySeed& seed)const
00612 {
00613   PTrajectoryStateOnDet pState( seed.startingState());
00614   const GeomDet* gdet  = (&(*tracker))->idToDet(DetId(pState.detId()));
00615   TSOS  State= trajectoryStateTransform::transientState( pState, &(gdet->surface()), 
00616                                            &(*magfield));
00617   return State;
00618 
00619 }
00620 
00621 void CRackTrajectoryBuilder::AddHit(Trajectory &traj,
00622                                      vector<const TrackingRecHit*>Hits, Propagator *currPropagator){
00623 
00624    if ( Hits.size() == 0 )
00625      return;
00626   
00627    if (debug_info) cout << "CRackTrajectoryBuilder::AddHit" << endl;
00628    if (debug_info) cout << "START " << traj.lastMeasurement().updatedState() << endl;
00629  
00630    vector <TrackingRecHitRange> hitRangeByDet;
00631    TrackingRecHitIterator prevDet;
00632  
00633    prevDet = Hits.begin();
00634    for( TrackingRecHitIterator iHit = Hits.begin(); iHit != Hits.end(); iHit++ )
00635      {
00636        if ( (*prevDet)->geographicalId() == (*iHit)->geographicalId() ) 
00637         continue;
00638         
00639        hitRangeByDet.push_back( make_pair( prevDet, iHit ) );
00640        prevDet = iHit;
00641       }
00642    hitRangeByDet.push_back( make_pair( prevDet, Hits.end() ) );
00643  
00645  
00646      if (fastPropagation) {
00647        for( TrackingRecHitRangeIterator iHitRange = hitRangeByDet.begin(); iHitRange != hitRangeByDet.end(); iHitRange++ )
00648          {
00649             const TrackingRecHit *currHit = *(iHitRange->first);
00650           DetId      currDet = currHit->geographicalId();
00651   
00652            TSOS prSt= currPropagator->propagate(traj.lastMeasurement().updatedState(),
00653                                               tracker->idToDet(currDet)->surface());
00654   
00655            if ( !prSt.isValid())
00656             {
00657               if (debug_info) cout << "Not Valid: PRST" << prSt.globalPosition();
00658  //           if (debug_info) cout << "Not Valid: HIT" << *currHit;
00659  
00660  
00661               continue;
00662             }
00663  
00664            TransientTrackingRecHit::RecHitPointer bestHit = RHBuilder->build( currHit );
00665            double chi2min = theEstimator->estimate( prSt, *bestHit).second;
00666  
00667            if (debug_info) cout << "Size " << iHitRange->first - (*iHitRange).second << endl;
00668            for( TrackingRecHitIterator iHit = (*iHitRange).first+1; iHit != iHitRange->second; iHit++ )
00669              {
00670                if (debug_info) cout << "loop3 " <<" "<< Hits.end() - iHit << endl;
00671  
00672                TransientTrackingRecHit::RecHitPointer tmpHit = RHBuilder->build( *iHit );
00673                double currChi2 = theEstimator->estimate(prSt, *tmpHit).second;
00674                if ( currChi2 < chi2min )
00675                  {
00676                    currChi2 = chi2min;
00677                    bestHit = tmpHit;
00678                  }
00679              }
00680            //now we have check if the track can be added to the trajectory
00681            if (debug_info) cout << chi2min << endl;
00682            if (chi2min < chi2cut)
00683              {
00684                if (debug_info) cout << "chi2 fine : " << chi2min << endl;
00685                TSOS UpdatedState= theUpdator->update( prSt, *bestHit );
00686                if (UpdatedState.isValid()){
00687                 hits.push_back(&(*bestHit));
00688                  traj.push( TM(prSt,UpdatedState, bestHit, chi2min) );
00689                 if (debug_info) edm::LogInfo("CosmicTrackFinder") <<
00690                    "STATE UPDATED WITH HIT AT POSITION "
00691                   
00692                                               << bestHit->globalPosition()
00693                                                <<UpdatedState<<" "
00694                                                <<traj.chiSquared();
00695                 if (debug_info) cout  <<
00696                    "STATE UPDATED WITH HIT AT POSITION "
00697                   
00698                                               << bestHit->globalPosition()
00699                                                <<UpdatedState<<" "
00700                                                <<traj.chiSquared();
00701                  if (debug_info) cout << "State is valid ..." << endl;
00702                 break; // now we need to 
00703                }
00704                else
00705                  {
00706                   edm::LogWarning("CosmicTrackFinder")<<" State can not be updated with hit at position " << endl;
00707                   TSOS UpdatedState= theUpdator->update( prSt, *bestHit );
00708                   if (UpdatedState.isValid()){
00709                     cout  <<
00710                       "NOT! UPDATED WITH HIT AT POSITION "
00711                       
00712                           << bestHit->globalPosition()
00713                           <<UpdatedState<<" "
00714                           <<traj.chiSquared();
00715                     
00716                   }
00717                 }
00718             } 
00719          }
00720      } //simple version end
00721      else 
00722      {
00723        //first sort the dets in the order they are traversed by the trajectory
00724        // we need three loops:
00725        // 1: loop as long as there can be an new hit added to the trajectory
00726        // 2: loop over all dets that might be hit
00727        // 3: loop over all hits on a certain det
00728        
00729        
00730        std::vector < std::pair<TrackingRecHitRangeIterator, TSOS> > trackHitCandidates;
00731        std::vector <std::pair<TrackingRecHitRangeIterator, TSOS> >::iterator iHitRange;
00732        std::vector <uint32_t> processedDets;
00733        do
00734         {
00735           
00736           //create vector of possibly hit detectors...
00737           trackHitCandidates.clear();      
00738           DetId currDet;
00739           for( TrackingRecHitRangeIterator iHit = hitRangeByDet.begin(); iHit != hitRangeByDet.end(); iHit++ )
00740             {
00741               const TrackingRecHit *currHit = *(iHit->first);
00742               currDet = currHit->geographicalId();
00743  
00744               if ( find(processedDets.begin(), processedDets.end(), currDet.rawId()) != processedDets.end() )
00745                 continue;
00746               
00747               TSOS prSt= currPropagator->propagate(traj.lastMeasurement().updatedState(),
00748                                                 tracker->idToDet(currDet)->surface());
00749               if ( ( !prSt.isValid() ) ||  (theEstimator->Chi2MeasurementEstimatorBase::estimate(prSt,tracker->idToDet(currDet)->surface() ) == false) )
00750  //           if ( ( !prSt.isValid() ) ||  (theEstimator->estimate(prSt,tracker->idToDet(currDet)->surface() ) == false) )  
00751                 continue;
00752      
00753               trackHitCandidates.push_back(  make_pair(iHit, prSt) );
00754             }
00755           
00756           if (!trackHitCandidates.size())
00757           break;
00758  
00759           if (debug_info) cout << Hits.size() << " (int) trackHitCandidates.begin() " << trackHitCandidates.size() << endl;
00760           if (debug_info) cout << "Before sorting ... " << endl; 
00761  
00762           if (debug_info)         
00763           for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ )
00764             {
00765               if (debug_info) cout << (tracker->idToDet((*(iHitRange->first->first))->geographicalId()))->position();
00766             }
00767           if (debug_info) cout << endl;
00768           
00769  
00770           stable_sort( trackHitCandidates.begin(), trackHitCandidates.end(), CompareDetByTraj(traj.lastMeasurement().updatedState())  );
00771           
00772           if (debug_info) cout << "After sorting ... " << endl; 
00773           if (debug_info)
00774             {
00775             for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ )
00776               {
00777                 if (debug_info) cout << (tracker->idToDet((*(iHitRange->first->first))->geographicalId()))->position();
00778               }
00779             cout << endl;
00780             }
00781  
00782         for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ )      //loop over dets
00783           {
00784             
00785             //now find the best hit of the detector
00786             if (debug_info) cout << "loop2 " << trackHitCandidates.size()  <<" " << trackHitCandidates.end() - iHitRange << endl;
00787             const TrackingRecHit *currHit = *(iHitRange->first->first);
00788             
00789             TransientTrackingRecHit::RecHitPointer bestHit = RHBuilder->build( currHit );
00790             TSOS currPrSt = (*iHitRange).second;
00791         
00792             if (debug_info) cout << "curr position" << bestHit->globalPosition();
00793             for( TrackingRecHitIterator iHit = (*iHitRange).first->first+1; iHit != iHitRange->first->second; iHit++ )
00794               {
00795                 TransientTrackingRecHit::RecHitPointer tmpHit = RHBuilder->build( *iHit );
00796                 if (debug_info) cout << "curr position" << tmpHit->globalPosition() ;
00797           
00798               }
00799           }
00800         if (debug_info) cout << "Cross check end ..." << endl;
00801  
00802  
00803         //just a simple test if the same hit can be added twice ...
00804         //      for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ )      //loop over all hits
00805         
00806         //      break;
00807         
00808         for( iHitRange = trackHitCandidates.begin(); iHitRange != trackHitCandidates.end(); iHitRange++ )      //loop over detsall hits
00809           {
00810             
00811             //now find the best hit of the detector
00812             if (debug_info) cout << "loop2 " << trackHitCandidates.size()  <<" " << trackHitCandidates.end() - iHitRange << endl;
00813             
00814             const TrackingRecHit *currHit = *(iHitRange->first->first);
00815  
00816             processedDets.push_back(currHit->geographicalId().rawId());
00817             
00818             
00819             TransientTrackingRecHit::RecHitPointer bestHit = RHBuilder->build( currHit );
00820  
00821             if (debug_info) cout << "curr position A" << bestHit->globalPosition() << endl;
00822             TSOS currPrSt = (*iHitRange).second;
00823             double chi2min = theEstimator->estimate( currPrSt, *bestHit).second;
00824  
00825             if (debug_info) cout << "Size " << iHitRange->first->second - (*iHitRange).first->first << endl; 
00826             for( TrackingRecHitIterator iHit = (*iHitRange).first->first+1; iHit != iHitRange->first->second; iHit++ )
00827               {
00828                 if (debug_info) cout << "loop3 " <<" "<< Hits.end() - iHit << endl;
00829               
00830                 TransientTrackingRecHit::RecHitPointer tmpHit = RHBuilder->build( *iHit );
00831                 if (debug_info) cout << "curr position B" << tmpHit->globalPosition() << endl;
00832                 double currChi2 = theEstimator->estimate(currPrSt, *tmpHit).second;
00833                 if ( currChi2 < chi2min )
00834                   {
00835                     if (debug_info) cout << "Is best hit" << endl;
00836                     chi2min = currChi2;
00837                     bestHit = tmpHit;
00838                   }
00839               }
00840             //now we have checked the det and can remove the entry from the vector...
00841           
00842             //if (debug_info) cout << "before erase ..." << endl;
00843             //this is to slow ...
00844             // hitRangeByDet.erase( (*iHitRange).first,(*iHitRange).first+1 );    
00845             //if (debug_info) cout << "after erase ..." << endl;
00846  
00847             if (debug_info) cout << chi2min << endl;
00848             //if the track can be added to the trajectory 
00849             if (chi2min < chi2cut)
00850               {
00851                 if (debug_info) cout << "chi2 fine : " << chi2min << endl;
00852  
00853                 //  if (debug_info) cout << "previaous state " << traj.lastMeasurement().updatedState() <<endl;
00854                 TSOS UpdatedState= theUpdator->update( currPrSt, *bestHit );
00855                 if (UpdatedState.isValid()){
00856  
00857                   hits.push_back(&(*bestHit));
00858                   traj.push( TM(currPrSt,UpdatedState, bestHit, chi2min) );
00859                   if (debug_info) edm::LogInfo("CosmicTrackFinder") <<
00860                     "STATE UPDATED WITH HIT AT POSITION "
00861                     //   <<tmphitbestdet->globalPosition()
00862                                                 <<UpdatedState<<" "
00863                                                 <<traj.chiSquared();
00864                   if (debug_info) cout << "Added Hit" << bestHit->globalPosition() << endl;
00865                   if (debug_info) cout << "State is valid ..." << UpdatedState << endl;
00866                   //cout << "updated state " << traj.lastMeasurement().updatedState() <<endl;
00867  
00868                   //          return; //break;
00869                   //
00870                   //          TSOS prSt= currPropagator->propagate(traj.lastMeasurement().updatedState(),
00871                   //                                          tracker->idToDet( bestHit->geographicalId() )->surface());
00872                   //      
00873                   //          if ( prSt.isValid())
00874                   //            cout << "the same hit can be added twice ..." << endl;  
00875                   //
00876                   
00877                   break;
00878                 }
00879                 else
00880                   {
00881                     if (debug_info) edm::LogWarning("CosmicTrackFinder")<<" State can not be updated with hit at position "
00882                                                         << bestHit->globalPosition();  
00883                     //          cout << "State can not be updated with hit at " << bestHit->globalPosition() << endl;
00884                   }
00885                 continue;
00886               }
00887             else
00888               {
00889                 //      cout << "chi2 to big : " << chi2min << endl;
00890               }
00891             if (debug_info) cout << " continue 1 " << endl;
00892           }
00893         //now we remove all already processed dets from the list ...
00894         // hitRangeByDet.erase( (*trackHitCandidates.begin()).first,(*iHitRange).first+1 );       
00895  
00896         if (debug_info) cout << " continue 2 " << endl;
00897        }
00898      //if this was the last exit
00899      while ( iHitRange != trackHitCandidates.end() );
00900      }
00901  
00902  
00903 }
00904 
00905 
00906 bool 
00907 CRackTrajectoryBuilder::qualityFilter(Trajectory traj){
00908   int ngoodhits=0;
00909   if(geometry=="MTCC"){
00910     std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> > hits= traj.recHits();
00911     std::vector< ConstReferenceCountingPointer< TransientTrackingRecHit> >::const_iterator hit;
00912     for(hit=hits.begin();hit!=hits.end();hit++){
00913       unsigned int iid=(*hit)->hit()->geographicalId().rawId();
00914       //CHECK FOR 3 hits r-phi
00915       if(((iid>>0)&0x3)!=1) ngoodhits++;
00916     }
00917   }
00918   else ngoodhits=traj.foundHits();
00919   
00920   if ( ngoodhits >= theMinHits) {
00921     return true;
00922   }
00923   else {
00924     return false;
00925   }
00926 }
00927 
00928 //----------------------------------------------------------------------------------------------------------
00929 // little helper function that returns false if both hits conatin the same information
00930 
00931 bool CRackTrajectoryBuilder::isDifferentStripReHit2D  (const SiStripRecHit2D& hitA, const  SiStripRecHit2D& hitB )
00932 {
00933   if ( hitA.geographicalId() != hitB.geographicalId() )
00934     return true;
00935   if ( hitA.localPosition().x() != hitB.localPosition().x() )
00936     return true;
00937   if ( hitA.localPosition().y() != hitB.localPosition().y() )
00938     return true;
00939   if ( hitA.localPosition().z() != hitB.localPosition().z() )
00940     return true;
00941 
00942   //  if (debug_info) cout << hitA.localPosition() << endl;
00943   //  if (debug_info) cout << hitB << endl;
00944 
00945   return false;
00946 }
00947 
00948 
00949 
00950 //----backfitting
00951 std::pair<TrajectoryStateOnSurface, const GeomDet*> 
00952 CRackTrajectoryBuilder::innerState( const Trajectory& traj) const
00953 {
00954   int lastFitted = 999;
00955   int nhits = traj.foundHits();
00956   if (nhits < lastFitted+1) lastFitted = nhits-1;
00957 
00958   std::vector<TrajectoryMeasurement> measvec = traj.measurements();
00959   TransientTrackingRecHit::ConstRecHitContainer firstHits;
00960 
00961   bool foundLast = false;
00962   int actualLast = -99;
00963   for (int i=lastFitted; i >= 0; i--) {
00964     if(measvec[i].recHit()->isValid()){
00965       if(!foundLast){
00966         actualLast = i; 
00967         foundLast = true;
00968       }
00969       firstHits.push_back( measvec[i].recHit());
00970     }
00971   }
00972   TSOS unscaledState = measvec[actualLast].updatedState();
00973   AlgebraicSymMatrix55 C=ROOT::Math::SMatrixIdentity();
00974   // C *= 100.;
00975 
00976   TSOS startingState( unscaledState.localParameters(), LocalTrajectoryError(C),
00977                       unscaledState.surface(),
00978                       thePropagator->magneticField());
00979 
00980   // cout << endl << "FitTester starts with state " << startingState << endl;
00981 
00982   KFTrajectoryFitter backFitter( *thePropagator,
00983                                  KFUpdator(),
00984                                  Chi2MeasurementEstimator( 100., 3));
00985 
00986   PropagationDirection backFitDirection = traj.direction() == alongMomentum ? oppositeToMomentum: alongMomentum;
00987 
00988   // only direction matters in this contest
00989   TrajectorySeed fakeSeed = TrajectorySeed(PTrajectoryStateOnDet() , 
00990                                            edm::OwnVector<TrackingRecHit>(),
00991                                            backFitDirection);
00992 
00993   vector<Trajectory> fitres = backFitter.fit( fakeSeed, firstHits, startingState);
00994 
00995   if (fitres.size() != 1) {
00996     // cout << "FitTester: first hits fit failed!" << endl;
00997     return std::pair<TrajectoryStateOnSurface, const GeomDet*>();
00998   }
00999 
01000   TrajectoryMeasurement firstMeas = fitres[0].lastMeasurement();
01001   TSOS firstState = firstMeas.updatedState();
01002 
01003   //  cout << "FitTester: Fitted first state " << firstState << endl;
01004   //cout << "FitTester: chi2 = " << fitres[0].chiSquared() << endl;
01005 
01006   TSOS initialState( firstState.localParameters(), LocalTrajectoryError(C),
01007                      firstState.surface(),
01008                      thePropagator->magneticField());
01009 
01010   return std::pair<TrajectoryStateOnSurface, const GeomDet*>( initialState, 
01011                                                               firstMeas.recHit()->det());
01012 }