CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoTracker/TrackProducer/src/TrackProducerAlgorithm.cc

Go to the documentation of this file.
00001 #include "RecoTracker/TrackProducer/interface/TrackProducerAlgorithm.h"
00002 
00003 #include "DataFormats/Common/interface/OrphanHandle.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 
00006 #include "MagneticField/Engine/interface/MagneticField.h"
00007 #include "Geometry/CommonDetUnit/interface/TrackingGeometry.h"
00008 
00009 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
00010 #include "DataFormats/TrackCandidate/interface/TrackCandidate.h"
00011 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
00012 
00013 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00014 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00015 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00016 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00017 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
00018 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00019 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00020 #include "RecoTracker/TrackProducer/interface/TrackingRecHitLessFromGlobalPosition.h"
00021 
00022 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00023 #include "FWCore/Utilities/interface/Exception.h"
00024 
00025 #include "RecoTracker/TransientTrackingRecHit/interface/TRecHit2DPosConstraint.h"
00026 #include "RecoTracker/TransientTrackingRecHit/interface/TRecHit1DMomConstraint.h"
00027 // #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00028 #include "DataFormats/GeometryCommonDetAlgo/interface/ErrorFrameTransformer.h"
00029 #include "TrackingTools/TrackFitters/interface/RecHitSorter.h"
00030 #include "DataFormats/TrackReco/interface/TrackBase.h"
00031 
00032 namespace {
00033 #ifdef STAT_TSB
00034   struct StatCount {
00035     long long totTrack=0;
00036     long long totLoop=0;
00037     long long totGsfTrack=0;
00038     long long totFound=0;
00039     long long totLost=0;
00040     long long totAlgo[12];
00041     void track(int l) {
00042       if (l>0) ++totLoop; else ++totTrack;
00043     }
00044     void hits(int f, int l) { totFound+=f; totLost+=l;} 
00045     void gsf() {++totGsfTrack;}
00046     void algo(int a) { if (a>=0 && a<12) ++totAlgo[a];}
00047 
00048 
00049     void print() const {
00050       std::cout << "TrackProducer stat\nTrack/Loop/Gsf/FoundHits/LostHits/algos "
00051                 <<  totTrack <<'/'<< totLoop <<'/'<< totGsfTrack  <<'/'<< totFound  <<'/'<< totLost;
00052       for (auto a : totAlgo) std::cout << '/'<< a;
00053         std::cout  << std::endl;
00054     }
00055     StatCount() {}
00056     ~StatCount() { print();}
00057   };
00058 
00059 #else
00060   struct StatCount {
00061     void track(int){}
00062     void hits(int, int){}
00063     void gsf(){}
00064     void algo(int){}
00065   };
00066 #endif
00067 
00068   StatCount statCount;
00069 
00070 }
00071 
00072 
00073 
00074 
00075 template <> bool
00076 TrackProducerAlgorithm<reco::Track>::buildTrack (const TrajectoryFitter * theFitter,
00077                                                  const Propagator * thePropagator,
00078                                                  AlgoProductCollection& algoResults,
00079                                                  TransientTrackingRecHit::RecHitContainer& hits,
00080                                                  TrajectoryStateOnSurface& theTSOS,
00081                                                  const TrajectorySeed& seed,
00082                                                  float ndof,
00083                                                  const reco::BeamSpot& bs,
00084                                                  SeedRef seedRef,
00085                                                  int qualityMask,signed char nLoops)                                             
00086 {
00087   //variable declarations
00088   reco::Track * theTrack;
00089   Trajectory * theTraj; 
00090   PropagationDirection seedDir = seed.direction();
00091       
00092   //perform the fit: the result's size is 1 if it succeded, 0 if fails
00093   Trajectory && trajTmp = theFitter->fitOne(seed, hits, theTSOS,(nLoops>0) ? TrajectoryFitter::looper : TrajectoryFitter::standard);
00094   if unlikely(!trajTmp.isValid()) return false;
00095   
00096   
00097   
00098   theTraj = new Trajectory(std::move(trajTmp));
00099   theTraj->setSeedRef(seedRef);
00100   
00101   statCount.hits(theTraj->foundHits(),theTraj->lostHits());
00102   statCount.algo(int(algo_));
00103 
00104   // TrajectoryStateOnSurface innertsos;
00105   // if (theTraj->direction() == alongMomentum) {
00106   //  innertsos = theTraj->firstMeasurement().updatedState();
00107   // } else { 
00108   //  innertsos = theTraj->lastMeasurement().updatedState();
00109   // }
00110   
00111   ndof = 0;
00112   for (auto const & tm : theTraj->measurements()) {
00113     auto const & h = tm.recHitR();
00114     if (h.isValid()) ndof = ndof + float(h.dimension())*h.weight();  // two virtual calls!
00115   }
00116   
00117   ndof -= 5.f;
00118   if unlikely(std::abs(theTSOS.magneticField()->nominalValue())<DBL_MIN) ++ndof;  // same as -4
00119  
00120  
00121   //if geometricInnerState_ is false the state for projection to beam line is the state attached to the first hit: to be used for loopers
00122   //if geometricInnerState_ is true the state for projection to beam line is the one from the (geometrically) closest measurement to the beam line: to be sued for non-collision tracks
00123   //the two shouuld give the same result for collision tracks that are NOT loopers
00124   TrajectoryStateOnSurface stateForProjectionToBeamLineOnSurface;
00125   if (geometricInnerState_) {
00126     stateForProjectionToBeamLineOnSurface = theTraj->closestMeasurement(GlobalPoint(bs.x0(),bs.y0(),bs.z0())).updatedState();
00127   } else {
00128     if (theTraj->direction() == alongMomentum) {
00129       stateForProjectionToBeamLineOnSurface = theTraj->firstMeasurement().updatedState();
00130     } else { 
00131       stateForProjectionToBeamLineOnSurface = theTraj->lastMeasurement().updatedState();
00132     }
00133   }
00134 
00135   if unlikely(!stateForProjectionToBeamLineOnSurface.isValid()){
00136     edm::LogError("CannotPropagateToBeamLine")<<"the state on the closest measurement isnot valid. skipping track.";
00137     delete theTraj;
00138     return false;
00139   }
00140   const FreeTrajectoryState & stateForProjectionToBeamLine=*stateForProjectionToBeamLineOnSurface.freeState();
00141   
00142   LogDebug("TrackProducer") << "stateForProjectionToBeamLine=" << stateForProjectionToBeamLine;
00143   
00144   TSCBLBuilderNoMaterial tscblBuilder;
00145   TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(stateForProjectionToBeamLine,bs);
00146   
00147   if unlikely(!tscbl.isValid()) {
00148     delete theTraj;
00149     return false;
00150   }
00151   
00152   GlobalPoint v = tscbl.trackStateAtPCA().position();
00153   math::XYZPoint  pos( v.x(), v.y(), v.z() );
00154   GlobalVector p = tscbl.trackStateAtPCA().momentum();
00155   math::XYZVector mom( p.x(), p.y(), p.z() );
00156   
00157   LogDebug("TrackProducer") << "pos=" << v << " mom=" << p << " pt=" << p.perp() << " mag=" << p.mag();
00158   
00159   theTrack = new reco::Track(theTraj->chiSquared(),
00160                              int(ndof),//FIXME fix weight() in TrackingRecHit
00161                              pos, mom, tscbl.trackStateAtPCA().charge(), 
00162                              tscbl.trackStateAtPCA().curvilinearError(),
00163                              algo_);
00164   
00165   theTrack->setQualityMask(qualityMask);
00166   theTrack->setNLoops(nLoops);
00167   
00168   LogDebug("TrackProducer") << "theTrack->pt()=" << theTrack->pt();
00169   
00170   LogDebug("TrackProducer") <<"track done\n";
00171   
00172   AlgoProduct aProduct(theTraj,std::make_pair(theTrack,seedDir));
00173   algoResults.push_back(aProduct);
00174   
00175   statCount.track(nLoops);
00176 
00177   return true;
00178 } 
00179 
00180 template <> bool
00181 TrackProducerAlgorithm<reco::GsfTrack>::buildTrack (const TrajectoryFitter * theFitter,
00182                                                     const Propagator * thePropagator,
00183                                                     AlgoProductCollection& algoResults,
00184                                                     TransientTrackingRecHit::RecHitContainer& hits,
00185                                                     TrajectoryStateOnSurface& theTSOS,
00186                                                     const TrajectorySeed& seed,
00187                                                     float ndof,
00188                                                     const reco::BeamSpot& bs,
00189                                                     SeedRef seedRef,
00190                                                     int qualityMask,signed char nLoops)
00191 {
00192   //variable declarations
00193   reco::GsfTrack * theTrack;
00194   Trajectory * theTraj; 
00195   PropagationDirection seedDir = seed.direction();
00196   
00197   Trajectory && trajTmp = theFitter->fitOne(seed, hits, theTSOS,(nLoops>0) ? TrajectoryFitter::looper: TrajectoryFitter::standard);
00198   if unlikely(!trajTmp.isValid()) return false;
00199   
00200   
00201   theTraj = new Trajectory( std::move(trajTmp) );
00202   theTraj->setSeedRef(seedRef);
00203   
00204   //  TrajectoryStateOnSurface innertsos;
00205   // TrajectoryStateOnSurface outertsos;
00206 
00207   // if (theTraj->direction() == alongMomentum) {
00208   //  innertsos = theTraj->firstMeasurement().updatedState();
00209   //  outertsos = theTraj->lastMeasurement().updatedState();
00210   // } else { 
00211   //  innertsos = theTraj->lastMeasurement().updatedState();
00212   //  outertsos = theTraj->firstMeasurement().updatedState();
00213   // }
00214   //     std::cout
00215   //       << "Nr. of first / last states = "
00216   //       << innertsos.components().size() << " "
00217   //       << outertsos.components().size() << std::endl;
00218   //     std::vector<TrajectoryStateOnSurface> components = 
00219   //       innertsos.components();
00220   //     double sinTheta = 
00221   //       sin(innertsos.globalMomentum().theta());
00222   //     for ( std::vector<TrajectoryStateOnSurface>::const_iterator ic=components.begin();
00223   //      ic!=components.end(); ic++ ) {
00224   //       std::cout << " comp " << ic-components.begin() << " "
00225   //            << (*ic).weight() << " "
00226   //            << (*ic).localParameters().vector()[0]/sinTheta << " "
00227   //            << sqrt((*ic).localError().matrix()[0][0])/sinTheta << std::endl;
00228   //     }
00229   
00230   ndof = 0;
00231   for (auto const & tm : theTraj->measurements()) {
00232     auto const & h = tm.recHitR();
00233     if (h.isValid()) ndof = ndof + h.dimension()*h.weight();
00234   }
00235   
00236   ndof = ndof - 5;
00237   if unlikely(std::abs(theTSOS.magneticField()->nominalValue())<DBL_MIN) ++ndof;  // same as -4
00238   
00239   
00240   //if geometricInnerState_ is false the state for projection to beam line is the state attached to the first hit: to be used for loopers
00241   //if geometricInnerState_ is true the state for projection to beam line is the one from the (geometrically) closest measurement to the beam line: to be sued for non-collision tracks
00242   //the two shouuld give the same result for collision tracks that are NOT loopers
00243   TrajectoryStateOnSurface stateForProjectionToBeamLineOnSurface;
00244   if (geometricInnerState_) {
00245     stateForProjectionToBeamLineOnSurface = theTraj->closestMeasurement(GlobalPoint(bs.x0(),bs.y0(),bs.z0())).updatedState();
00246   } else {
00247     if (theTraj->direction() == alongMomentum) {
00248       stateForProjectionToBeamLineOnSurface = theTraj->firstMeasurement().updatedState();
00249     } else { 
00250       stateForProjectionToBeamLineOnSurface = theTraj->lastMeasurement().updatedState();
00251     }
00252   }
00253 
00254   if unlikely(!stateForProjectionToBeamLineOnSurface.isValid()){
00255       edm::LogError("CannotPropagateToBeamLine")<<"the state on the closest measurement isnot valid. skipping track.";
00256       delete theTraj;
00257       return false;
00258     }    
00259   
00260   const FreeTrajectoryState & stateForProjectionToBeamLine=*stateForProjectionToBeamLineOnSurface.freeState();
00261   
00262   LogDebug("GsfTrackProducer") << "stateForProjectionToBeamLine=" << stateForProjectionToBeamLine;
00263   
00264   TSCBLBuilderNoMaterial tscblBuilder;
00265   TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(stateForProjectionToBeamLine,bs);
00266   
00267   if unlikely(tscbl.isValid()==false) {
00268       delete theTraj;
00269       return false;
00270     }
00271   
00272   GlobalPoint v = tscbl.trackStateAtPCA().position();
00273   math::XYZPoint  pos( v.x(), v.y(), v.z() );
00274   GlobalVector p = tscbl.trackStateAtPCA().momentum();
00275   math::XYZVector mom( p.x(), p.y(), p.z() );
00276   
00277   LogDebug("GsfTrackProducer") << "pos=" << v << " mom=" << p << " pt=" << p.perp() << " mag=" << p.mag();
00278   
00279   theTrack = new reco::GsfTrack(theTraj->chiSquared(),
00280                                 int(ndof),//FIXME fix weight() in TrackingRecHit
00281                                 //                             theTraj->foundHits(),//FIXME to be fixed in Trajectory.h
00282                                 //                             0, //FIXME no corresponding method in trajectory.h
00283                                 //                             theTraj->lostHits(),//FIXME to be fixed in Trajectory.h
00284                                 pos, mom, tscbl.trackStateAtPCA().charge(), tscbl.trackStateAtPCA().curvilinearError());    
00285   theTrack->setAlgorithm(algo_);
00286   
00287   LogDebug("GsfTrackProducer") <<"track done\n";
00288   
00289   AlgoProduct aProduct(theTraj,std::make_pair(theTrack,seedDir));
00290   LogDebug("GsfTrackProducer") <<"track done1\n";
00291   algoResults.push_back(aProduct);
00292   LogDebug("GsfTrackProducer") <<"track done2\n";
00293   
00294   statCount.gsf();
00295   return true;
00296 }