CMS 3D CMS Logo

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