00001 #include "RecoMuon/StandAloneTrackFinder/interface/ExhaustiveMuonTrajectoryBuilder.h" 00002 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00003 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h" 00004 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHitBuilder.h" 00005 00006 00007 ExhaustiveMuonTrajectoryBuilder::ExhaustiveMuonTrajectoryBuilder(const edm::ParameterSet & pset, 00008 const MuonServiceProxy* proxy) 00009 : theTrajBuilder(pset, proxy), 00010 theSeeder(), 00011 theService(proxy) 00012 { 00013 } 00014 00015 00016 ExhaustiveMuonTrajectoryBuilder::~ExhaustiveMuonTrajectoryBuilder() {} 00017 00018 00019 MuonTrajectoryBuilder::TrajectoryContainer 00020 ExhaustiveMuonTrajectoryBuilder::trajectories(const TrajectorySeed& seed) 00021 { 00022 LocalTrajectoryParameters localTrajectoryParameters(seed.startingState().parameters()); 00023 LocalVector p(localTrajectoryParameters.momentum()); 00024 int rawId = seed.startingState().detId(); 00025 DetId detId(rawId); 00026 bool isBarrel = (detId.subdetId() == 1); 00027 // homemade local-to-global 00028 double pt = (isBarrel) ? -p.z() : p.perp(); 00029 pt *= localTrajectoryParameters.charge(); 00030 float err00 = seed.startingState().error(0); 00031 // float p_err = sqr(sptmean/(ptmean*ptmean)); 00032 // mat[0][0]= p_err; 00033 float sigmapt = sqrt(err00)*pt*pt; 00034 TrajectorySeed::range range = seed.recHits(); 00035 TrajectoryContainer result; 00036 // Make a new seed based on each segment, using the original pt and sigmapt 00037 for(TrajectorySeed::const_iterator recHitItr = range.first; 00038 recHitItr != range.second; ++recHitItr) 00039 { 00040 const GeomDet * geomDet = theService->trackingGeometry()->idToDet((*recHitItr).geographicalId()); 00041 MuonTransientTrackingRecHit::MuonRecHitPointer muonRecHit 00042 = MuonTransientTrackingRecHit::specificBuild(geomDet,&*recHitItr); 00043 TrajectorySeed tmpSeed(theSeeder.createSeed(pt, sigmapt, muonRecHit)); 00044 TrajectoryContainer trajectories(theTrajBuilder.trajectories(tmpSeed)); 00045 result.insert(result.end(), trajectories.begin(), trajectories.end()); 00046 } 00047 // choose the best trajectory 00048 if(!result.empty()) clean(result); 00049 return result; 00050 } 00051 00052 00053 MuonTrajectoryBuilder::CandidateContainer 00054 ExhaustiveMuonTrajectoryBuilder::trajectories(const TrackCand&) 00055 { 00056 return CandidateContainer(); 00057 } 00058 00059 00060 00061 void ExhaustiveMuonTrajectoryBuilder::setEvent(const edm::Event& event) 00062 { 00063 theTrajBuilder.setEvent(event); 00064 } 00065 00066 00067 void ExhaustiveMuonTrajectoryBuilder::clean(TrajectoryContainer & trajectories) const 00068 { 00069 // choose the one with the most hits, and the smallest chi-squared 00070 int best_nhits = 0; 00071 unsigned best = 0; 00072 unsigned ntraj = trajectories.size(); 00073 for(unsigned i = 0; i < ntraj; ++i) 00074 { 00075 int nhits = trajectories[i]->foundHits(); 00076 if(nhits > best_nhits) 00077 { 00078 best_nhits = nhits; 00079 best = i; 00080 } 00081 else if(nhits == best_nhits && trajectories[i]->chiSquared() < trajectories[best]->chiSquared()) 00082 { 00083 best = i; 00084 } 00085 } 00086 TrajectoryContainer result; 00087 for(unsigned i = 0; i < ntraj; ++i) 00088 { 00089 if(i == best) { 00090 result.push_back(trajectories[best]); 00091 } else { 00092 delete trajectories[i]; 00093 } 00094 } 00095 trajectories.swap(result); 00096 } 00097