Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:14

0001 #include "RecoMuon/StandAloneTrackFinder/interface/ExhaustiveMuonTrajectoryBuilder.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
0004 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHitBuilder.h"
0005 
0006 ExhaustiveMuonTrajectoryBuilder::ExhaustiveMuonTrajectoryBuilder(const edm::ParameterSet& pset,
0007                                                                  const MuonServiceProxy* proxy,
0008                                                                  edm::ConsumesCollector& iC)
0009     : theTrajBuilder(pset, proxy, iC), theSeeder(), theService(proxy) {}
0010 
0011 ExhaustiveMuonTrajectoryBuilder::~ExhaustiveMuonTrajectoryBuilder() {}
0012 
0013 MuonTrajectoryBuilder::TrajectoryContainer ExhaustiveMuonTrajectoryBuilder::trajectories(const TrajectorySeed& seed) {
0014   LocalTrajectoryParameters localTrajectoryParameters(seed.startingState().parameters());
0015   LocalVector p(localTrajectoryParameters.momentum());
0016   int rawId = seed.startingState().detId();
0017   DetId detId(rawId);
0018   bool isBarrel = (detId.subdetId() == 1);
0019   // homemade local-to-global
0020   double pt = (isBarrel) ? -p.z() : p.perp();
0021   pt *= localTrajectoryParameters.charge();
0022   float err00 = seed.startingState().error(0);
0023   //   float p_err = sqr(sptmean/(ptmean*ptmean));
0024   //  mat[0][0]= p_err;
0025   float sigmapt = sqrt(err00) * pt * pt;
0026   TrajectoryContainer result;
0027   // Make a new seed based on each segment, using the original pt and sigmapt
0028   for (auto const& recHit : seed.recHits()) {
0029     const GeomDet* geomDet = theService->trackingGeometry()->idToDet(recHit.geographicalId());
0030     auto muonRecHit = MuonTransientTrackingRecHit::specificBuild(geomDet, &recHit);
0031     TrajectorySeed tmpSeed(theSeeder.createSeed(pt, sigmapt, muonRecHit));
0032     TrajectoryContainer trajectories(theTrajBuilder.trajectories(tmpSeed));
0033     result.insert(
0034         result.end(), std::make_move_iterator(trajectories.begin()), std::make_move_iterator(trajectories.end()));
0035   }
0036   return result;
0037 }
0038 
0039 MuonTrajectoryBuilder::CandidateContainer ExhaustiveMuonTrajectoryBuilder::trajectories(const TrackCand&) {
0040   return CandidateContainer();
0041 }
0042 
0043 void ExhaustiveMuonTrajectoryBuilder::setEvent(const edm::Event& event) { theTrajBuilder.setEvent(event); }
0044 
0045 void ExhaustiveMuonTrajectoryBuilder::clean(TrajectoryContainer& trajectories) const {
0046   // choose the one with the most hits, and the smallest chi-squared
0047   if (trajectories.empty()) {
0048     return;
0049   }
0050   int best_nhits = 0;
0051   unsigned best = 0;
0052   unsigned ntraj = trajectories.size();
0053   for (unsigned i = 0; i < ntraj; ++i) {
0054     int nhits = trajectories[i]->foundHits();
0055     if (nhits > best_nhits) {
0056       best_nhits = nhits;
0057       best = i;
0058     } else if (nhits == best_nhits && trajectories[i]->chiSquared() < trajectories[best]->chiSquared()) {
0059       best = i;
0060     }
0061   }
0062   TrajectoryContainer result;
0063   result.reserve(1);
0064   result.emplace_back(std::move(trajectories[best]));
0065   trajectories.swap(result);
0066 }