CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoMuon/GlobalTrackingTools/src/GlobalTrajectoryBuilderBase.cc

Go to the documentation of this file.
00001 
00026 #include "RecoMuon/GlobalTrackingTools/interface/GlobalTrajectoryBuilderBase.h"
00027 
00028 //---------------
00029 // C++ Headers --
00030 //---------------
00031 
00032 #include <iostream>
00033 #include <algorithm>
00034 
00035 //-------------------------------
00036 // Collaborating Class Headers --
00037 //-------------------------------
00038 
00039 #include "FWCore/Framework/interface/Event.h"
00040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00041 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00042 
00043 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00044 #include "TrackingTools/TrackFitters/interface/RecHitLessByDet.h"
00045 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00046 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00047 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00048 #include "TrackingTools/TrackRefitter/interface/TrackTransformer.h"
00049 
00050 #include "DataFormats/Math/interface/deltaR.h"
00051 
00052 #include "DataFormats/DetId/interface/DetId.h"
00053 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
00054 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
00055 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
00056 
00057 #include "DataFormats/TrackReco/interface/Track.h"
00058 #include "DataFormats/TrackReco/interface/TrackExtraFwd.h"
00059 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00060 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00061 
00062 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00063 
00064 #include "RecoMuon/GlobalTrackingTools/interface/GlobalMuonTrackMatcher.h"
00065 #include "RecoMuon/GlobalTrackingTools/interface/GlobalMuonRefitter.h"
00066 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHitBuilder.h"
00067 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHit.h"
00068 #include "RecoMuon/TrackingTools/interface/MuonCandidate.h"
00069 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
00070 
00071 #include "RecoMuon/GlobalTrackingTools/interface/MuonTrackingRegionBuilder.h"
00072 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00073 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00074 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00075 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00076 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00077 
00078 #include "RecoTracker/TkTrackingRegions/interface/TkTrackingRegionsMargin.h"
00079 #include "RecoTracker/TkMSParametrization/interface/PixelRecoRange.h"
00080 
00081 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00082 #include "RecoTracker/TransientTrackingRecHit/interface/TSiStripRecHit2DLocalPos.h"
00083 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00084 
00085 using namespace std;
00086 using namespace edm;
00087 
00088 //----------------
00089 // Constructors --
00090 //----------------
00091 GlobalTrajectoryBuilderBase::GlobalTrajectoryBuilderBase(const edm::ParameterSet& par,
00092                                                          const MuonServiceProxy* service) : 
00093   theTrackMatcher(0),theLayerMeasurements(0),theTrackTransformer(0),theRegionBuilder(0), theService(service),theGlbRefitter(0) {
00094 
00095   theCategory = par.getUntrackedParameter<string>("Category", "Muon|RecoMuon|GlobalMuon|GlobalTrajectoryBuilderBase");
00096 
00097   
00098   ParameterSet trackMatcherPSet = par.getParameter<ParameterSet>("GlobalMuonTrackMatcher");
00099   theTrackMatcher = new GlobalMuonTrackMatcher(trackMatcherPSet,theService);
00100   
00101   theTrackerPropagatorName = par.getParameter<string>("TrackerPropagator");
00102 
00103   ParameterSet trackTransformerPSet = par.getParameter<ParameterSet>("TrackTransformer");
00104   theTrackTransformer = new TrackTransformer(trackTransformerPSet);
00105 
00106   ParameterSet regionBuilderPSet = par.getParameter<ParameterSet>("MuonTrackingRegionBuilder");
00107 
00108   theRegionBuilder = new MuonTrackingRegionBuilder(regionBuilderPSet,theService);
00109 
00110   // TrackRefitter parameters
00111   ParameterSet refitterParameters = par.getParameter<ParameterSet>("GlbRefitterParameters");
00112   theGlbRefitter = new GlobalMuonRefitter(refitterParameters, theService);
00113 
00114   theMuonHitsOption = refitterParameters.getParameter<int>("MuonHitsOption");
00115 
00116   theTrackerRecHitBuilderName = par.getParameter<string>("TrackerRecHitBuilder");
00117   theMuonRecHitBuilderName = par.getParameter<string>("MuonRecHitBuilder");  
00118 
00119   theRPCInTheFit = par.getParameter<bool>("RefitRPCHits");
00120 
00121   theTECxScale = par.getParameter<double>("ScaleTECxFactor");
00122   theTECyScale = par.getParameter<double>("ScaleTECyFactor");
00123   thePtCut = par.getParameter<double>("PtCut");
00124   thePCut = par.getParameter<double>("PCut");
00125 
00126   theCacheId_TRH = 0;
00127 
00128 }
00129 
00130 
00131 //--------------
00132 // Destructor --
00133 //--------------
00134 GlobalTrajectoryBuilderBase::~GlobalTrajectoryBuilderBase() {
00135 
00136   if (theTrackMatcher) delete theTrackMatcher;
00137   if (theRegionBuilder) delete theRegionBuilder;
00138   if (theTrackTransformer) delete theTrackTransformer;
00139   if (theGlbRefitter) delete theGlbRefitter;
00140 }
00141 
00142 
00143 //
00144 // set Event
00145 //
00146 void GlobalTrajectoryBuilderBase::setEvent(const edm::Event& event) {
00147   
00148   theEvent = &event;
00149 
00150   theTrackTransformer->setServices(theService->eventSetup());
00151   theRegionBuilder->setEvent(event);
00152 
00153   theGlbRefitter->setEvent(event);
00154   theGlbRefitter->setServices(theService->eventSetup());
00155 
00156   unsigned long long newCacheId_TRH = theService->eventSetup().get<TransientRecHitRecord>().cacheIdentifier();
00157   if ( newCacheId_TRH != theCacheId_TRH ) {
00158     LogDebug(theCategory) << "TransientRecHitRecord changed!";
00159     theCacheId_TRH = newCacheId_TRH;
00160     theService->eventSetup().get<TransientRecHitRecord>().get(theTrackerRecHitBuilderName,theTrackerRecHitBuilder);
00161     theService->eventSetup().get<TransientRecHitRecord>().get(theMuonRecHitBuilderName,theMuonRecHitBuilder);
00162   }
00163 
00164   //Retrieve tracker topology from geometry
00165   edm::ESHandle<TrackerTopology> tTopoHand;
00166   theService->eventSetup().get<IdealGeometryRecord>().get(tTopoHand);
00167   tTopo_=tTopoHand.product();
00168 
00169 }
00170 
00171 
00172 //
00173 // build a combined tracker-muon trajectory
00174 //
00175 MuonCandidate::CandidateContainer 
00176 GlobalTrajectoryBuilderBase::build(const TrackCand& staCand,
00177                                    MuonCandidate::CandidateContainer& tkTrajs ) const {
00178 
00179   LogTrace(theCategory) << " Begin Build" << endl;
00180 
00181   // tracker trajectory should be built and refit before this point
00182   if ( tkTrajs.empty() ) return CandidateContainer();
00183 
00184   // add muon hits and refit/smooth trajectories
00185   CandidateContainer refittedResult;
00186   ConstRecHitContainer muonRecHits = getTransientRecHits(*(staCand.second));
00187 
00188   // check order of muon measurements
00189   if ( (muonRecHits.size() > 1) &&
00190        ( muonRecHits.front()->globalPosition().mag() >
00191          muonRecHits.back()->globalPosition().mag() ) ) {
00192     LogTrace(theCategory)<< "   reverse order: ";
00193   }
00194 
00195   for ( CandidateContainer::const_iterator it = tkTrajs.begin(); it != tkTrajs.end(); it++ ) {
00196 
00197     // cut on tracks with low momenta
00198     LogTrace(theCategory)<< "   Track p and pT " << (*it)->trackerTrack()->p() << " " << (*it)->trackerTrack()->pt();
00199     if(  (*it)->trackerTrack()->p() < thePCut || (*it)->trackerTrack()->pt() < thePtCut  ) continue;
00200 
00201     ConstRecHitContainer trackerRecHits;
00202     if ((*it)->trackerTrack().isNonnull()) {
00203       trackerRecHits = getTransientRecHits(*(*it)->trackerTrack());
00204     } else {
00205       LogDebug(theCategory)<<"     NEED HITS FROM TRAJ";
00206       //trackerRecHits = (*it)->trackerTrajectory()->recHits();
00207     }
00208 
00209     // check for single TEC RecHits in trajectories in the overalp region
00210     if ( fabs((*it)->trackerTrack()->eta()) > 0.95 && fabs((*it)->trackerTrack()->eta()) < 1.15 && (*it)->trackerTrack()->pt() < 60 ) {
00211       if ( theTECxScale < 0 || theTECyScale < 0 )
00212         trackerRecHits = selectTrackerHits(trackerRecHits);
00213       else
00214         fixTEC(trackerRecHits,theTECxScale,theTECyScale);
00215     }
00216                   
00217     RefitDirection recHitDir = checkRecHitsOrdering(trackerRecHits);
00218     if ( recHitDir == outToIn ) reverse(trackerRecHits.begin(),trackerRecHits.end());
00219 
00220     reco::TransientTrack tTT((*it)->trackerTrack(),&*theService->magneticField(),theService->trackingGeometry());
00221     TrajectoryStateOnSurface innerTsos = tTT.innermostMeasurementState();
00222 
00223     edm::RefToBase<TrajectorySeed> tmpSeed;
00224     if((*it)->trackerTrack()->seedRef().isAvailable()) tmpSeed = (*it)->trackerTrack()->seedRef();
00225 
00226     if ( !innerTsos.isValid() ) {
00227       LogTrace(theCategory) << "     inner Trajectory State is invalid. ";
00228       continue;
00229     }
00230   
00231     innerTsos.rescaleError(100.);
00232                   
00233     TC refitted0,refitted1;
00234     MuonCandidate* finalTrajectory = 0;
00235     Trajectory *tkTrajectory = 0;
00236 
00237     // tracker only track
00238     if ( ! ((*it)->trackerTrajectory() && (*it)->trackerTrajectory()->isValid()) ) { 
00239       refitted0 = theTrackTransformer->transform((*it)->trackerTrack()) ;
00240       if (!refitted0.empty()) tkTrajectory = new Trajectory(*(refitted0.begin())); 
00241       else LogWarning(theCategory)<< "     Failed to load tracker track trajectory";
00242     } else tkTrajectory = (*it)->trackerTrajectory();
00243     if (tkTrajectory) tkTrajectory->setSeedRef(tmpSeed);
00244 
00245     // full track with all muon hits using theGlbRefitter    
00246     ConstRecHitContainer allRecHits = trackerRecHits;
00247     allRecHits.insert(allRecHits.end(), muonRecHits.begin(),muonRecHits.end());
00248     refitted1 = theGlbRefitter->refit( *(*it)->trackerTrack(), tTT, allRecHits,theMuonHitsOption, tTopo_);
00249     LogTrace(theCategory)<<"     This track-sta refitted to " << refitted1.size() << " trajectories";
00250 
00251     Trajectory *glbTrajectory1 = 0;
00252     if (!refitted1.empty()) glbTrajectory1 = new Trajectory(*(refitted1.begin()));
00253     else LogDebug(theCategory)<< "     Failed to load global track trajectory 1"; 
00254     if (glbTrajectory1) glbTrajectory1->setSeedRef(tmpSeed);
00255     
00256     finalTrajectory = 0;
00257     if(glbTrajectory1 && tkTrajectory) finalTrajectory = new MuonCandidate(glbTrajectory1, (*it)->muonTrack(), (*it)->trackerTrack(), 
00258                                         tkTrajectory? new Trajectory(*tkTrajectory) : 0);
00259 
00260     if ( finalTrajectory ) 
00261       refittedResult.push_back(finalTrajectory);
00262      
00263     if(tkTrajectory) delete tkTrajectory;
00264   }
00265 
00266   // choose the best global fit for this Standalone Muon based on the track probability
00267   CandidateContainer selectedResult;
00268   MuonCandidate* tmpCand = 0;
00269   if ( refittedResult.size() > 0 ) tmpCand = *(refittedResult.begin());
00270   double minProb = 9999;
00271 
00272   for (CandidateContainer::const_iterator iter=refittedResult.begin(); iter != refittedResult.end(); iter++) {
00273     double prob = trackProbability(*(*iter)->trajectory());
00274     LogTrace(theCategory)<<"   refitted-track-sta with pT " << (*iter)->trackerTrack()->pt() << " has probability " << prob;
00275 
00276     if (prob < minProb) {
00277       minProb = prob;
00278       tmpCand = (*iter);
00279     }
00280   }
00281 
00282   if ( tmpCand )  selectedResult.push_back(new MuonCandidate(new Trajectory(*(tmpCand->trajectory())), tmpCand->muonTrack(), tmpCand->trackerTrack(), 
00283                                                              (tmpCand->trackerTrajectory())? new Trajectory( *(tmpCand->trackerTrajectory()) ):0 ) );
00284 
00285   for (CandidateContainer::const_iterator it = refittedResult.begin(); it != refittedResult.end(); ++it) {
00286     if ( (*it)->trajectory() ) delete (*it)->trajectory();
00287     if ( (*it)->trackerTrajectory() ) delete (*it)->trackerTrajectory();
00288     if ( *it ) delete (*it);
00289   }
00290   refittedResult.clear();
00291 
00292   return selectedResult;
00293 
00294 }
00295 
00296 
00297 //
00298 // select tracker tracks within a region of interest
00299 //
00300 vector<GlobalTrajectoryBuilderBase::TrackCand> 
00301 GlobalTrajectoryBuilderBase::chooseRegionalTrackerTracks(const TrackCand& staCand, 
00302                                                          const vector<TrackCand>& tkTs) {
00303   
00304   // define eta-phi region
00305   RectangularEtaPhiTrackingRegion regionOfInterest = defineRegionOfInterest(staCand.second);
00306   
00307   // get region's etaRange and phiMargin
00308   //UNUSED:  PixelRecoRange<float> etaRange = regionOfInterest.etaRange();
00309   //UNUSED:  TkTrackingRegionsMargin<float> phiMargin = regionOfInterest.phiMargin();
00310 
00311   vector<TrackCand> result;
00312 
00313   double deltaR_max = 1.0;
00314 
00315   for ( vector<TrackCand>::const_iterator is = tkTs.begin(); is != tkTs.end(); ++is ) {
00316     // check if each trackCand is in region of interest
00317 //    bool inEtaRange = etaRange.inside(is->second->eta());
00318 //    bool inPhiRange = (fabs(Geom::Phi<float>(is->second->phi()) - Geom::Phi<float>(regionOfInterest.direction().phi())) < phiMargin.right() ) ? true : false ;
00319 
00320     double deltaR_tmp = deltaR(static_cast<double>(regionOfInterest.direction().eta()),
00321                                                            static_cast<double>(regionOfInterest.direction().phi()),
00322                                                            is->second->eta(), is->second->phi());
00323 
00324     // for each trackCand in region, add trajectory and add to result
00325     //if ( inEtaRange && inPhiRange ) {
00326     if (deltaR_tmp < deltaR_max) {
00327       TrackCand tmpCand = TrackCand(*is);
00328       result.push_back(tmpCand);
00329     }
00330   }
00331 
00332   return result; 
00333 
00334 }
00335 
00336 
00337 //
00338 // define a region of interest within the tracker
00339 //
00340 RectangularEtaPhiTrackingRegion 
00341 GlobalTrajectoryBuilderBase::defineRegionOfInterest(const reco::TrackRef& staTrack) const {
00342 
00343   RectangularEtaPhiTrackingRegion* region1 = theRegionBuilder->region(staTrack);
00344   
00345   TkTrackingRegionsMargin<float> etaMargin(fabs(region1->etaRange().min() - region1->etaRange().mean()),
00346                                            fabs(region1->etaRange().max() - region1->etaRange().mean()));
00347   
00348   RectangularEtaPhiTrackingRegion region2(region1->direction(),
00349                                           region1->origin(),
00350                                           region1->ptMin(),
00351                                           region1->originRBound(),
00352                                           region1->originZBound(),
00353                                           etaMargin,
00354                                           region1->phiMargin());
00355   
00356   delete region1;
00357   return region2;
00358   
00359 }
00360 
00361 
00362 //
00363 // calculate the tail probability (-ln(P)) of a fit
00364 //
00365 double 
00366 GlobalTrajectoryBuilderBase::trackProbability(const Trajectory& track) const {
00367 
00368   if ( track.ndof() > 0 && track.chiSquared() > 0 ) { 
00369     return -LnChiSquaredProbability(track.chiSquared(), track.ndof());
00370   } else {
00371     return 0.0;
00372   }
00373 
00374 }
00375 
00376 
00377 //
00378 // print RecHits
00379 //
00380 void GlobalTrajectoryBuilderBase::printHits(const ConstRecHitContainer& hits) const {
00381 
00382   LogTrace(theCategory) << "Used RecHits: " << hits.size();
00383   for (ConstRecHitContainer::const_iterator ir = hits.begin(); ir != hits.end(); ir++ ) {
00384     if ( !(*ir)->isValid() ) {
00385       LogTrace(theCategory) << "invalid RecHit";
00386       continue; 
00387     }
00388     
00389     const GlobalPoint& pos = (*ir)->globalPosition();
00390     
00391     LogTrace(theCategory) 
00392       << "r = " << sqrt(pos.x() * pos.x() + pos.y() * pos.y())
00393       << "  z = " << pos.z()
00394       << "  dimension = " << (*ir)->dimension()
00395       << "  " << (*ir)->det()->geographicalId().det()
00396       << "  " << (*ir)->det()->subDetector();
00397 
00398   }
00399 
00400 }
00401 
00402 //
00403 // check order of RechIts on a trajectory
00404 //
00405 GlobalTrajectoryBuilderBase::RefitDirection
00406 GlobalTrajectoryBuilderBase::checkRecHitsOrdering(const TransientTrackingRecHit::ConstRecHitContainer& recHits) const {
00407 
00408   if ( !recHits.empty() ) {
00409     ConstRecHitContainer::const_iterator frontHit = recHits.begin();
00410     ConstRecHitContainer::const_iterator backHit  = recHits.end() - 1;
00411     while ( !(*frontHit)->isValid() && frontHit != backHit ) {frontHit++;}
00412     while ( !(*backHit)->isValid() && backHit != frontHit )  {backHit--;}
00413 
00414     double rFirst = (*frontHit)->globalPosition().mag();
00415     double rLast  = (*backHit) ->globalPosition().mag();
00416 
00417     if ( rFirst < rLast ) return inToOut;
00418     else if (rFirst > rLast) return outToIn;
00419     else {
00420       LogError(theCategory) << "Impossible to determine the rechits order" << endl;
00421       return undetermined;
00422     }
00423   }
00424   else {
00425     LogError(theCategory) << "Impossible to determine the rechits order" << endl;
00426     return undetermined;
00427   }
00428 }
00429 
00430 
00431 //
00432 // select trajectories with only a single TEC hit
00433 //
00434 GlobalTrajectoryBuilderBase::ConstRecHitContainer 
00435 GlobalTrajectoryBuilderBase::selectTrackerHits(const ConstRecHitContainer& all) const {
00436  
00437   int nTEC(0);
00438 
00439   ConstRecHitContainer hits;
00440   for (ConstRecHitContainer::const_iterator i = all.begin(); i != all.end(); i++) {
00441     if ( !(*i)->isValid() ) continue;
00442     if ( (*i)->det()->geographicalId().det() == DetId::Tracker &&
00443          (*i)->det()->geographicalId().subdetId() == StripSubdetector::TEC) {
00444       nTEC++;
00445     } else {
00446       hits.push_back((*i).get());
00447     }
00448     if ( nTEC > 1 ) return all;
00449   }
00450 
00451   return hits;
00452 
00453 }
00454 
00455 
00456 //
00457 // rescale errors of outermost TEC RecHit
00458 //
00459 void GlobalTrajectoryBuilderBase::fixTEC(ConstRecHitContainer& all,
00460                                          double scl_x,
00461                                          double scl_y) const {
00462 
00463   int nTEC(0);
00464   ConstRecHitContainer::iterator lone_tec;
00465 
00466   for ( ConstRecHitContainer::iterator i = all.begin(); i != all.end(); i++) {
00467     if ( !(*i)->isValid() ) continue;
00468     
00469     if ( (*i)->det()->geographicalId().det() == DetId::Tracker &&
00470          (*i)->det()->geographicalId().subdetId() == StripSubdetector::TEC) {
00471       lone_tec = i;
00472       nTEC++;
00473       
00474       if ( (i+1) != all.end() && (*(i+1))->isValid() &&
00475           (*(i+1))->det()->geographicalId().det() == DetId::Tracker &&
00476           (*(i+1))->det()->geographicalId().subdetId() == StripSubdetector::TEC) {
00477         nTEC++;
00478         break;
00479       }
00480     }
00481     
00482     if (nTEC > 1) break;
00483   }
00484   
00485   int hitDet = (*lone_tec)->hit()->geographicalId().det();
00486   int hitSubDet = (*lone_tec)->hit()->geographicalId().subdetId();
00487   if ( nTEC == 1 && (*lone_tec)->hit()->isValid() &&
00488        hitDet == DetId::Tracker && hitSubDet == StripSubdetector::TEC) {
00489     
00490     // rescale the TEC rechit error matrix in its rotated frame
00491     const SiStripRecHit2D* strip = dynamic_cast<const SiStripRecHit2D*>((*lone_tec)->hit());
00492     const TSiStripRecHit2DLocalPos* Tstrip = dynamic_cast<const TSiStripRecHit2DLocalPos*>((*lone_tec).get());
00493     if (strip && Tstrip->det() && Tstrip) {
00494       LocalPoint pos = Tstrip->localPosition();
00495       if ((*lone_tec)->detUnit()) {
00496         const StripTopology* topology = dynamic_cast<const StripTopology*>(&(*lone_tec)->detUnit()->topology());
00497         if (topology) {
00498           // rescale the local error along/perp the strip by a factor
00499           float angle = topology->stripAngle(topology->strip((*lone_tec)->hit()->localPosition()));
00500           LocalError error = Tstrip->localPositionError();
00501           LocalError rotError = error.rotate(angle);
00502           LocalError scaledError(rotError.xx() * scl_x * scl_x, 0, rotError.yy() * scl_y * scl_y);
00503           error = scaledError.rotate(-angle);
00504           MuonTransientTrackingRecHit* mtt_rechit;
00505           if (strip->cluster().isNonnull()) {
00509             SiStripRecHit2D* st = new SiStripRecHit2D(pos,error,
00510                                                       (*lone_tec)->geographicalId().rawId(),
00511                                                       strip->cluster());
00512             *lone_tec = mtt_rechit->build((*lone_tec)->det(),st);
00513           }
00514           else {
00515             SiStripRecHit2D* st = new SiStripRecHit2D(pos,error,
00516                                                       (*lone_tec)->geographicalId().rawId(),
00517                                                       strip->cluster_regional());
00518             *lone_tec = mtt_rechit->build((*lone_tec)->det(),st);
00519           }
00520         }
00521       }
00522     }
00523   }
00524 
00525 }
00526 
00527 
00528 //
00529 // get transient RecHits
00530 //
00531 TransientTrackingRecHit::ConstRecHitContainer
00532 GlobalTrajectoryBuilderBase::getTransientRecHits(const reco::Track& track) const {
00533 
00534   TransientTrackingRecHit::ConstRecHitContainer result;
00535   
00536   
00537 
00538   TrajectoryStateOnSurface currTsos = trajectoryStateTransform::innerStateOnSurface(track, *theService->trackingGeometry(), &*theService->magneticField());
00539 
00540   for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
00541     if((*hit)->isValid()) {
00542       DetId recoid = (*hit)->geographicalId();
00543       if ( recoid.det() == DetId::Tracker ) {
00544         TransientTrackingRecHit::RecHitPointer ttrhit = theTrackerRecHitBuilder->build(&**hit);
00545         if (!ttrhit->hit()->hasPositionAndError()){
00546           TrajectoryStateOnSurface predTsos =  theService->propagator(theTrackerPropagatorName)->propagate(currTsos, theService->trackingGeometry()->idToDet(recoid)->surface());
00547           
00548           if ( !predTsos.isValid() ) {
00549             edm::LogError("MissingTransientHit")
00550               <<"Could not get a tsos on the hit surface. We will miss a tracking hit.";
00551             continue; 
00552           }
00553           currTsos = predTsos;
00554           TransientTrackingRecHit::RecHitPointer preciseHit = ttrhit->clone(predTsos);
00555           result.push_back(preciseHit);
00556         }else{
00557           result.push_back(ttrhit);
00558         }
00559       } else if ( recoid.det() == DetId::Muon ) {
00560         if ( (*hit)->geographicalId().subdetId() == 3 && !theRPCInTheFit) {
00561           LogDebug(theCategory) << "RPC Rec Hit discarded"; 
00562           continue;
00563         }
00564         result.push_back(theMuonRecHitBuilder->build(&**hit));
00565       }
00566     }
00567   }
00568   
00569   return result;
00570 }