Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002  * \class TrackerSeedCleaner
0003  *  Reference class for seeds cleaning
0004  *  Seeds Cleaner based on sharedHits cleaning, direction cleaning and pt cleaning
0005     \author A. Grelli -  Purdue University, Pavia University
0006  */
0007 
0008 #include "RecoMuon/TrackerSeedGenerator/interface/TrackerSeedCleaner.h"
0009 
0010 //---------------
0011 // C++ Headers --
0012 //---------------
0013 #include <vector>
0014 
0015 //-------------------------------
0016 // Collaborating Class Headers --
0017 //-------------------------------
0018 #include "FWCore/Framework/interface/Event.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 
0021 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
0022 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0023 #include "DataFormats/Math/interface/deltaPhi.h"
0024 
0025 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
0026 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
0027 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
0028 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
0029 #include "RecoTracker/TkTrackingRegions/interface/TkTrackingRegionsMargin.h"
0030 #include "RecoTracker/TkMSParametrization/interface/PixelRecoRange.h"
0031 
0032 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
0033 
0034 using namespace std;
0035 using namespace edm;
0036 
0037 //
0038 // inizialization
0039 //
0040 void TrackerSeedCleaner::init(const MuonServiceProxy* service) {
0041   theProxyService = service;
0042 
0043   theRedundantCleaner = new RedundantSeedCleaner();
0044 }
0045 
0046 //
0047 //
0048 //
0049 void TrackerSeedCleaner::setEvent(const edm::Event& event) { event.getByToken(beamspotToken_, bsHandle_); }
0050 
0051 //
0052 // clean seeds
0053 //
0054 void TrackerSeedCleaner::clean(const reco::TrackRef& muR,
0055                                const RectangularEtaPhiTrackingRegion& region,
0056                                tkSeeds& seeds) {
0057   // call the shared input cleaner
0058   if (cleanBySharedHits)
0059     theRedundantCleaner->define(seeds);
0060 
0061   theTTRHBuilder = theProxyService->eventSetup().getHandle(theTTRHBuilderToken);
0062 
0063   LogDebug("TrackerSeedCleaner") << seeds.size() << " trajectory seeds to the events before cleaning" << endl;
0064 
0065   //check the validity otherwise vertexing
0066   const reco::BeamSpot& bs = *bsHandle_;
0067   /*reco track and seeds as arguments. Seeds eta and phi are checked and 
0068    based on deviation from L2 eta and phi seed is accepted or not*/
0069 
0070   std::vector<TrajectorySeed> result;
0071 
0072   TSCBLBuilderNoMaterial tscblBuilder;
0073   // PerigeeConversions tspConverter;
0074   for (TrajectorySeedCollection::iterator seed = seeds.begin(); seed < seeds.end(); ++seed) {
0075     if (seed->nHits() < 2)
0076       continue;
0077     //get parameters and errors from the seed state
0078     TransientTrackingRecHit::RecHitPointer recHit = theTTRHBuilder->build(&*(seed->recHits().end() - 1));
0079     TrajectoryStateOnSurface state = trajectoryStateTransform::transientState(
0080         seed->startingState(), recHit->surface(), theProxyService->magneticField().product());
0081 
0082     TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed =
0083         tscblBuilder(*state.freeState(), bs);  //as in TrackProducerAlgorithms
0084     if (!tsAtClosestApproachSeed.isValid())
0085       continue;
0086     GlobalPoint vSeed1 = tsAtClosestApproachSeed.trackStateAtPCA().position();
0087     GlobalVector pSeed = tsAtClosestApproachSeed.trackStateAtPCA().momentum();
0088     GlobalPoint vSeed(vSeed1.x() - bs.x0(), vSeed1.y() - bs.y0(), vSeed1.z() - bs.z0());
0089 
0090     //eta,phi info from seeds
0091     double etaSeed = state.globalMomentum().eta();
0092     double phiSeed = pSeed.phi();
0093 
0094     //if the limits are too stringent rescale limits
0095     typedef PixelRecoRange<float> Range;
0096     typedef TkTrackingRegionsMargin<float> Margin;
0097 
0098     Range etaRange = region.etaRange();
0099     double etaLimit = (fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) < 0.1)
0100                           ? 0.1
0101                           : fabs(fabs(etaRange.max()) - fabs(etaRange.mean()));
0102 
0103     Margin phiMargin = region.phiMargin();
0104     double phiLimit = (phiMargin.right() < 0.1) ? 0.1 : phiMargin.right();
0105 
0106     double ptSeed = pSeed.perp();
0107     double ptMin = (region.ptMin() > 3.5) ? 3.5 : region.ptMin();
0108     // Clean
0109     bool inEtaRange = etaSeed >= (etaRange.mean() - etaLimit) && etaSeed <= (etaRange.mean() + etaLimit);
0110     bool inPhiRange = (fabs(deltaPhi(phiSeed, double(region.direction().phi()))) < phiLimit);
0111     // pt cleaner
0112     bool inPtRange = ptSeed >= ptMin && ptSeed <= 2 * (muR->pt());
0113 
0114     // save efficiency don't clean triplets with pt cleaner
0115     if (seed->nHits() == 3)
0116       inPtRange = true;
0117 
0118     // use pt and angle cleaners
0119     if (inPtRange && usePt_Cleaner && !useDirection_Cleaner) {
0120       result.push_back(*seed);
0121       LogDebug("TrackerSeedCleaner") << " Keeping the seed : this seed passed pt selection";
0122     }
0123 
0124     // use only angle default option
0125     if (inEtaRange && inPhiRange && !usePt_Cleaner && useDirection_Cleaner) {
0126       result.push_back(*seed);
0127       LogDebug("TrackerSeedCleaner") << " Keeping the seed : this seed passed direction selection";
0128     }
0129 
0130     // use all the cleaners
0131     if (inEtaRange && inPhiRange && inPtRange && usePt_Cleaner && useDirection_Cleaner) {
0132       result.push_back(*seed);
0133       LogDebug("TrackerSeedCleaner") << " Keeping the seed : this seed passed direction and pt selection";
0134     }
0135 
0136     LogDebug("TrackerSeedCleaner") << " eta for current seed " << etaSeed << "\n"
0137                                    << " phi for current seed " << phiSeed << "\n"
0138                                    << " eta for L2 track  " << muR->eta() << "\n"
0139                                    << " phi for L2 track  " << muR->phi() << "\n";
0140   }
0141 
0142   //the new seeds collection
0143   if (!result.empty() && (useDirection_Cleaner || usePt_Cleaner))
0144     seeds.swap(result);
0145 
0146   LogDebug("TrackerSeedCleaner") << seeds.size() << " trajectory seeds to the events after cleaning" << endl;
0147 
0148   return;
0149 }