CMS 3D CMS Logo

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/PatternTools/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/TSCPBuilderNoMaterial.h"
00023 #include "TrackingTools/PatternTools/interface/TrajectoryStateClosestToBeamLineBuilder.h"
00024 #include "Utilities/General/interface/CMSexception.h"
00025 
00026 #include "RecoTracker/TransientTrackingRecHit/interface/TRecHit2DPosConstraint.h"
00027 #include "RecoTracker/TransientTrackingRecHit/interface/TRecHit1DMomConstraint.h"
00028 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00029 #include "DataFormats/GeometryCommonDetAlgo/interface/ErrorFrameTransformer.h"
00030 #include "TrackingTools/TrackFitters/interface/RecHitSorter.h"
00031 #include "DataFormats/TrackReco/interface/TrackBase.h"
00032 
00033 template <> bool
00034 TrackProducerAlgorithm<reco::Track>::buildTrack (const TrajectoryFitter * theFitter,
00035                                                  const Propagator * thePropagator,
00036                                                  AlgoProductCollection& algoResults,
00037                                                  TransientTrackingRecHit::RecHitContainer& hits,
00038                                                  TrajectoryStateOnSurface& theTSOS,
00039                                                  const TrajectorySeed& seed,
00040                                                  float ndof,
00041                                                  const reco::BeamSpot& bs,
00042                                                  SeedRef seedRef)
00043 {
00044   //variable declarations
00045   std::vector<Trajectory> trajVec;
00046   reco::Track * theTrack;
00047   Trajectory * theTraj; 
00048   PropagationDirection seedDir = seed.direction();
00049       
00050   //perform the fit: the result's size is 1 if it succeded, 0 if fails
00051   trajVec = theFitter->fit(seed, hits, theTSOS);
00052   
00053   LogDebug("TrackProducer") <<" FITTER FOUND "<< trajVec.size() << " TRAJECTORIES" <<"\n";
00054   TrajectoryStateOnSurface innertsos;
00055   
00056   if (trajVec.size() != 0){
00057 
00058     theTraj = new Trajectory( trajVec.front() );
00059     theTraj->setSeedRef(seedRef);
00060     
00061     if (theTraj->direction() == alongMomentum) {
00062       innertsos = theTraj->firstMeasurement().updatedState();
00063     } else { 
00064       innertsos = theTraj->lastMeasurement().updatedState();
00065     }
00066     
00067     ndof = 0;
00068     TransientTrackingRecHit::RecHitContainer validHits;
00069     theTraj->validRecHits(validHits);
00070     for (TransientTrackingRecHit::RecHitContainer::iterator h=validHits.begin();h!=validHits.end();++h)
00071       ndof = ndof + ((*h)->dimension())*((*h)->weight());
00072     if (theTSOS.magneticField()->inTesla(GlobalPoint(0,0,0)).mag2()<DBL_MIN) ndof = ndof - 4;
00073     else ndof = ndof - 5;
00074     
00075     TrajectoryStateClosestToBeamLineBuilder tscblBuilder;
00076     //    const FreeTrajectoryState & stateForProjectionToBeamLine=*innertsos.freeState();
00077     const FreeTrajectoryState & stateForProjectionToBeamLine=*theTraj->closestMeasurement(GlobalPoint(bs.x0(),bs.y0(),bs.z0())).updatedState().freeState();
00078 
00079     LogDebug("TrackProducer") << "stateForProjectionToBeamLine=" << stateForProjectionToBeamLine;
00080 
00081     TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(stateForProjectionToBeamLine,bs);
00082 
00083     if (tscbl.isValid()==false) {
00084         delete theTraj;
00085         return false;
00086     }
00087 
00088     GlobalPoint v = tscbl.trackStateAtPCA().position();
00089     math::XYZPoint  pos( v.x(), v.y(), v.z() );
00090     GlobalVector p = tscbl.trackStateAtPCA().momentum();
00091     math::XYZVector mom( p.x(), p.y(), p.z() );
00092 
00093     LogDebug("TrackProducer") << "pos=" << v << " mom=" << p << " pt=" << p.perp() << " mag=" << p.mag();
00094 
00095     algoName_ = conf_.getParameter<std::string>( "AlgorithmName" );  
00096     reco::TrackBase::TrackAlgorithm Algo = reco::TrackBase::algoByName(algoName_);
00097     
00098     theTrack = new reco::Track(theTraj->chiSquared(),
00099                                int(ndof),//FIXME fix weight() in TrackingRecHit
00100                                pos, mom, tscbl.trackStateAtPCA().charge(), 
00101                                tscbl.trackStateAtPCA().curvilinearError(),
00102                                Algo);
00103     
00104     LogDebug("TrackProducer") << "theTrack->pt()=" << theTrack->pt();
00105 
00106     LogDebug("TrackProducer") <<"track done\n";
00107 
00108     AlgoProduct aProduct(theTraj,std::make_pair(theTrack,seedDir));
00109     algoResults.push_back(aProduct);
00110     
00111     return true;
00112   } 
00113   else  return false;
00114 }
00115 
00116 template <> bool
00117 TrackProducerAlgorithm<reco::GsfTrack>::buildTrack (const TrajectoryFitter * theFitter,
00118                                                     const Propagator * thePropagator,
00119                                                     AlgoProductCollection& algoResults,
00120                                                     TransientTrackingRecHit::RecHitContainer& hits,
00121                                                     TrajectoryStateOnSurface& theTSOS,
00122                                                     const TrajectorySeed& seed,
00123                                                     float ndof,
00124                                                     const reco::BeamSpot& bs,
00125                                                     SeedRef seedRef)
00126 {
00127   //variable declarations
00128   std::vector<Trajectory> trajVec;
00129   reco::GsfTrack * theTrack;
00130   Trajectory * theTraj; 
00131   PropagationDirection seedDir = seed.direction();
00132       
00133   //perform the fit: the result's size is 1 if it succeded, 0 if fails
00134   trajVec = theFitter->fit(seed, hits, theTSOS);
00135   
00136   LogDebug("GsfTrackProducer") <<" FITTER FOUND "<< trajVec.size() << " TRAJECTORIES" <<"\n";
00137   
00138   TrajectoryStateOnSurface innertsos;
00139   TrajectoryStateOnSurface outertsos;
00140   
00141   if (trajVec.size() != 0){
00142 
00143     theTraj = new Trajectory( trajVec.front() );
00144     theTraj->setSeedRef(seedRef);
00145     
00146     if (theTraj->direction() == alongMomentum) {
00147       innertsos = theTraj->firstMeasurement().updatedState();
00148       outertsos = theTraj->lastMeasurement().updatedState();
00149     } else { 
00150       innertsos = theTraj->lastMeasurement().updatedState();
00151       outertsos = theTraj->firstMeasurement().updatedState();
00152     }
00153 //     std::cout
00154 //       << "Nr. of first / last states = "
00155 //       << innertsos.components().size() << " "
00156 //       << outertsos.components().size() << std::endl;
00157 //     std::vector<TrajectoryStateOnSurface> components = 
00158 //       innertsos.components();
00159 //     double sinTheta = 
00160 //       sin(innertsos.globalMomentum().theta());
00161 //     for ( std::vector<TrajectoryStateOnSurface>::const_iterator ic=components.begin();
00162 //        ic!=components.end(); ic++ ) {
00163 //       std::cout << " comp " << ic-components.begin() << " "
00164 //              << (*ic).weight() << " "
00165 //              << (*ic).localParameters().vector()[0]/sinTheta << " "
00166 //              << sqrt((*ic).localError().matrix()[0][0])/sinTheta << std::endl;
00167 //     }
00168     
00169     ndof = 0;
00170     TransientTrackingRecHit::RecHitContainer validHits;
00171     theTraj->validRecHits(validHits);
00172     for (TransientTrackingRecHit::RecHitContainer::iterator h=validHits.begin();h!=validHits.end();++h)
00173       ndof = ndof + ((*h)->dimension())*((*h)->weight());
00174     if (theTSOS.magneticField()->inTesla(GlobalPoint(0,0,0)).mag2()<DBL_MIN) ndof = ndof - 4;
00175     else ndof = ndof - 5;
00176    
00177     TrajectoryStateClosestToBeamLineBuilder tscblBuilder;
00178     //    const FreeTrajectoryState & stateForProjectionToBeamLine=*innertsos.freeState();
00179     const FreeTrajectoryState & stateForProjectionToBeamLine=*theTraj->closestMeasurement(GlobalPoint(bs.x0(),bs.y0(),bs.z0())).updatedState().freeState();
00180 
00181     LogDebug("TrackProducer") << "stateForProjectionToBeamLine=" << stateForProjectionToBeamLine;
00182 
00183     TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(stateForProjectionToBeamLine,bs);
00184 
00185     if (tscbl.isValid()==false) {
00186         delete theTraj;
00187         return false;
00188     }
00189 
00190     GlobalPoint v = tscbl.trackStateAtPCA().position();
00191     math::XYZPoint  pos( v.x(), v.y(), v.z() );
00192     GlobalVector p = tscbl.trackStateAtPCA().momentum();
00193     math::XYZVector mom( p.x(), p.y(), p.z() );
00194 
00195     LogDebug("TrackProducer") << "pos=" << v << " mom=" << p << " pt=" << p.perp() << " mag=" << p.mag();
00196 
00197     theTrack = new reco::GsfTrack(theTraj->chiSquared(),
00198                                   int(ndof),//FIXME fix weight() in TrackingRecHit
00199                                   //                           theTraj->foundHits(),//FIXME to be fixed in Trajectory.h
00200                                   //                           0, //FIXME no corresponding method in trajectory.h
00201                                   //                           theTraj->lostHits(),//FIXME to be fixed in Trajectory.h
00202                                   pos, mom, tscbl.trackStateAtPCA().charge(), tscbl.trackStateAtPCA().curvilinearError());    
00203 
00204     LogDebug("GsfTrackProducer") <<"track done\n";
00205 
00206     AlgoProduct aProduct(theTraj,std::make_pair(theTrack,seedDir));
00207     LogDebug("GsfTrackProducer") <<"track done1\n";
00208     algoResults.push_back(aProduct);
00209     LogDebug("GsfTrackProducer") <<"track done2\n";
00210     
00211     return true;
00212   } 
00213   else  return false;
00214 }

Generated on Tue Jun 9 17:46:00 2009 for CMSSW by  doxygen 1.5.4