Go to the documentation of this file.00001
00009 #include "HLTrigger/Muon/interface/HLTMuonL3PreFilter.h"
00010
00011 #include "DataFormats/Common/interface/Handle.h"
00012 #include "DataFormats/Common/interface/RefToBase.h"
00013
00014 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00015 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
00016
00017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00018 #include "DataFormats/TrackReco/interface/Track.h"
00019 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00020 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
00021 #include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
00022 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00023 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
00024 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
00025
00026
00027
00028
00029 using namespace std;
00030 using namespace edm;
00031 using namespace reco;
00032 using namespace trigger;
00033
00034 HLTMuonL3PreFilter::HLTMuonL3PreFilter(const ParameterSet& iConfig) :
00035 beamspotTag_ (iConfig.getParameter< edm::InputTag > ("BeamSpotTag")),
00036 candTag_ (iConfig.getParameter<InputTag > ("CandTag")),
00037 previousCandTag_ (iConfig.getParameter<InputTag > ("PreviousCandTag")),
00038 min_N_ (iConfig.getParameter<int> ("MinN")),
00039 max_Eta_ (iConfig.getParameter<double> ("MaxEta")),
00040 min_Nhits_ (iConfig.getParameter<int> ("MinNhits")),
00041 max_Dr_ (iConfig.getParameter<double> ("MaxDr")),
00042 max_Dz_ (iConfig.getParameter<double> ("MaxDz")),
00043 min_Pt_ (iConfig.getParameter<double> ("MinPt")),
00044 nsigma_Pt_ (iConfig.getParameter<double> ("NSigmaPt")),
00045 saveTag_ (iConfig.getUntrackedParameter<bool> ("SaveTag",false))
00046 {
00047
00048 LogDebug("HLTMuonL3PreFilter")
00049 << " CandTag/MinN/MaxEta/MinNhits/MaxDr/MaxDz/MinPt/NSigmaPt : "
00050 << candTag_.encode()
00051 << " " << min_N_
00052 << " " << max_Eta_
00053 << " " << min_Nhits_
00054 << " " << max_Dr_
00055 << " " << max_Dz_
00056 << " " << min_Pt_
00057 << " " << nsigma_Pt_;
00058
00059
00060 produces<TriggerFilterObjectWithRefs>();
00061 }
00062
00063 HLTMuonL3PreFilter::~HLTMuonL3PreFilter()
00064 {
00065 }
00066
00067
00068
00069
00070
00071
00072 bool
00073 HLTMuonL3PreFilter::filter(Event& iEvent, const EventSetup& iSetup)
00074 {
00075
00076
00077
00078
00079
00080
00081 auto_ptr<TriggerFilterObjectWithRefs>
00082 filterproduct (new TriggerFilterObjectWithRefs(path(),module()));
00083
00084
00085 Handle<RecoChargedCandidateCollection> mucands;
00086 if(saveTag_)filterproduct->addCollectionTag(candTag_);
00087 iEvent.getByLabel (candTag_,mucands);
00088
00089 std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;
00090 unsigned int maxI = mucands->size();
00091 for (unsigned int i=0;i!=maxI;++i){
00092 TrackRef tk = (*mucands)[i].track();
00093 edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef = tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
00094 TrackRef staTrack = l3seedRef->l2Track();
00095 LogDebug("HLTMuonL3PreFilter") <<"L2 from: "<<iEvent.getProvenance(staTrack.id()).moduleLabel() <<" index: "<<staTrack.key();
00096 L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands,i));
00097 }
00098
00099 Handle<TriggerFilterObjectWithRefs> previousLevelCands;
00100 iEvent.getByLabel (previousCandTag_,previousLevelCands);
00101 BeamSpot beamSpot;
00102 Handle<BeamSpot> recoBeamSpotHandle;
00103 iEvent.getByLabel(beamspotTag_,recoBeamSpotHandle);
00104 beamSpot = *recoBeamSpotHandle;
00105
00106
00107
00108 vector<RecoChargedCandidateRef> vl2cands;
00109 previousLevelCands->getObjects(TriggerMuon,vl2cands);
00110
00111
00112 int n = 0;
00113 std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > ::iterator L2toL3s_it = L2toL3s.begin();
00114 std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > ::iterator L2toL3s_end = L2toL3s.end();
00115 LogDebug("HLTMuonL3PreFilter")<<"looking at: "<<L2toL3s.size()<<" L2->L3s from: "<<mucands->size();
00116 for (; L2toL3s_it!=L2toL3s_end; ++L2toL3s_it){
00117
00118 if (!triggeredByLevel2(L2toL3s_it->first,vl2cands)) continue;
00119
00120
00121 unsigned int iTk=0;
00122 unsigned int maxItk=L2toL3s_it->second.size();
00123 for (; iTk!=maxItk; iTk++){
00124
00125 RecoChargedCandidateRef & cand=L2toL3s_it->second[iTk];
00126 TrackRef tk = cand->track();
00127
00128 LogDebug("HLTMuonL3PreFilter") << " Muon in loop, q*pt= " << tk->charge()*tk->pt() << ", eta= " << tk->eta() << ", hits= " << tk->numberOfValidHits() << ", d0= " << tk->d0() << ", dz= " << tk->dz();
00129
00130
00131 if (fabs(tk->eta())>max_Eta_) continue;
00132
00133
00134 if (tk->numberOfValidHits()<min_Nhits_) continue;
00135
00136
00137
00138 if (fabs(tk->dxy(beamSpot.position()))>max_Dr_) continue;
00139
00140
00141 if (fabs(tk->dz())>max_Dz_) continue;
00142
00143
00144 double pt = tk->pt();
00145 double err0 = tk->error(0);
00146 double abspar0 = fabs(tk->parameter(0));
00147 double ptLx = pt;
00148
00149 if (abspar0>0) ptLx += nsigma_Pt_*err0/abspar0*pt;
00150 LogTrace("HLTMuonL3PreFilter") << " ...Muon in loop, pt= "
00151 << pt << ", ptLx= " << ptLx;
00152 if (ptLx<min_Pt_) continue;
00153
00154 filterproduct->addObject(TriggerMuon,cand);
00155 n++;
00156 break;
00157 }
00158 }
00159
00160 vector<RecoChargedCandidateRef> vref;
00161 filterproduct->getObjects(TriggerMuon,vref);
00162 for (unsigned int i=0; i<vref.size(); i++ ) {
00163 RecoChargedCandidateRef candref = RecoChargedCandidateRef(vref[i]);
00164 TrackRef tk = candref->get<TrackRef>();
00165 LogDebug("HLTMuonL3PreFilter")
00166 << " Track passing filter: pt= " << tk->pt() << ", eta: "
00167 << tk->eta();
00168 }
00169
00170
00171 const bool accept (n >= min_N_);
00172
00173
00174 iEvent.put(filterproduct);
00175
00176 LogDebug("HLTMuonL3PreFilter") << " >>>>> Result of HLTMuonL3PreFilter is " << accept << ", number of muons passing thresholds= " << n;
00177
00178 return accept;
00179 }
00180 bool
00181 HLTMuonL3PreFilter::triggeredByLevel2(const TrackRef& staTrack,vector<RecoChargedCandidateRef>& vcands)
00182 {
00183 bool ok=false;
00184 for (unsigned int i=0; i<vcands.size(); i++) {
00185 if ( vcands[i]->get<TrackRef>() == staTrack ) {
00186 ok=true;
00187 LogDebug("HLTMuonL3PreFilter") << "The L2 track triggered";
00188 break;
00189 }
00190 }
00191 return ok;
00192 }
00193