CMS 3D CMS Logo

HLTMuonL3SimplePreFilter.cc
Go to the documentation of this file.
1 
9 
17 
19 
21 #include <iostream>
22 
23 //
24 // constructors and destructor
25 //
26  using namespace std;
27  using namespace edm;
28  using namespace trigger;
29  using namespace reco;
30 
32  candTag_ (iConfig.getParameter< edm::InputTag > ("CandTag") ),
33  previousCandTag_ (iConfig.getParameter< edm::InputTag > ("PreviousCandTag")),
34  beamspotTag_ (iConfig.getParameter< edm::InputTag > ("BeamSpotTag")),
35  min_N_ (iConfig.getParameter<int> ("MinN")),
36  max_Eta_ (iConfig.getParameter<double> ("MaxEta")),
37  min_Nhits_ (iConfig.getParameter<int> ("MinNhits")),
38  max_Dz_ (iConfig.getParameter<double> ("MaxDz")),
39  min_DxySig_ (iConfig.getParameter<double> ("MinDxySig")),
40  min_Pt_ (iConfig.getParameter<double> ("MinPt")),
41  nsigma_Pt_ (iConfig.getParameter<double> ("NSigmaPt")),
42  max_NormalizedChi2_ (iConfig.getParameter<double> ("MaxNormalizedChi2")),
43  max_DXYBeamSpot_ (iConfig.getParameter<double> ("MaxDXYBeamSpot")),
44  min_DXYBeamSpot_ (iConfig.getParameter<double> ("MinDXYBeamSpot")),
45  min_NmuonHits_ (iConfig.getParameter<int> ("MinNmuonHits")),
46  max_PtDifference_ (iConfig.getParameter<double> ("MaxPtDifference")),
47  min_TrackPt_ (iConfig.getParameter<double> ("MinTrackPt")),
48  matchPreviousCand_( iConfig.getParameter<bool>("MatchToPreviousCand") )
49 {
50  candToken_ = consumes<reco::RecoChargedCandidateCollection>(candTag_);
51  previousCandToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_);
52  beamspotToken_ = consumes<reco::BeamSpot>(beamspotTag_);
53 
54 }
55 
56 
58 
59 //
60 // member functions
61 //
62 void
66  desc.add<edm::InputTag>("BeamSpotTag",edm::InputTag("hltOfflineBeamSpot"));
67  desc.add<edm::InputTag>("CandTag",edm::InputTag("hltL3MuonCandidates"));
68  desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag(""));
69  desc.add<int>("MinN",1);
70  desc.add<double>("MaxEta",2.5);
71  desc.add<int>("MinNhits",0);
72  desc.add<double>("MaxDz",9999.0);
73  desc.add<double>("MinDxySig",-1.0);
74  desc.add<double>("MinPt",3.0);
75  desc.add<double>("NSigmaPt",0.0);
76  desc.add<double>("MaxNormalizedChi2",9999.0);
77  desc.add<double>("MaxDXYBeamSpot",9999.0);
78  desc.add<double>("MinDXYBeamSpot",-1.0);
79  desc.add<int>("MinNmuonHits",0);
80  desc.add<double>("MaxPtDifference",9999.0);
81  desc.add<double>("MinTrackPt",0.0);
82  desc.add<bool>("MatchToPreviousCand", true);
83  descriptions.add("hltMuonL3SimplePreFilter", desc);
84 }
85 
86 // ------------ method called to produce the data ------------
87 bool
89 {
90 
91  // All HLT filters must create and fill an HLT filter object,
92  // recording any reconstructed physics objects satisfying (or not)
93  // this HLT filter, and place it in the Event.
94 
95  if (saveTags()) filterproduct.addCollectionTag(candTag_);
96 
97  // get hold of trks
99  iEvent.getByToken (candToken_,mucands);
100  Handle<TriggerFilterObjectWithRefs> previousLevelCands;
101  iEvent.getByToken (previousCandToken_,previousLevelCands);
102  vector<RecoChargedCandidateRef> vcands;
103  if (previousLevelCands.isValid()) {
104  previousLevelCands->getObjects(TriggerMuon,vcands);
105  }
106 
107  Handle<BeamSpot> recoBeamSpotHandle;
108  iEvent.getByToken(beamspotToken_,recoBeamSpotHandle);
109 
110  // Number of objects passing the L3 Trigger:
111  int n = 0;
112  for (unsigned int iMu=0; iMu<mucands->size(); iMu++) {
113 
114  RecoChargedCandidateRef cand(mucands,iMu);
115  LogDebug("HLTMuonL3SimplePreFilter") << "cand isNonnull " << cand.isNonnull();
116 
117  //did this candidate triggered at previous stage.
118  if (matchPreviousCand_ && !triggerdByPreviousLevel(cand,vcands)) continue;
119 
120  if (std::abs(cand->eta())>max_Eta_) continue;
121 
122  TrackRef tk = cand->track();
123  LogDebug("HLTMuonL3SimplePreFilter") << " Muon in loop, q*pt= " << tk->charge()*tk->pt() <<" (" << cand->charge()*cand->pt() << ") "
124  << ", eta= " << tk->eta() << " (" << cand->eta() << ") " << ", hits= " << tk->numberOfValidHits()
125  << ", d0= " << tk->d0() << ", dz= " << tk->dz();
126 
127  // cut on number of hits
128  if (tk->numberOfValidHits()<min_Nhits_) continue;
129 
130  //normalizedChi2 cut
131  if (tk->normalizedChi2() > max_NormalizedChi2_ ) continue;
132 
133 
134  if (recoBeamSpotHandle.isValid()){
135  const BeamSpot& beamSpot = *recoBeamSpotHandle;
136 
137  //dz cut
138  if (std::abs((cand->vz()-beamSpot.z0()) - ((cand->vx()-beamSpot.x0())*cand->px()+(cand->vy()-beamSpot.y0())*cand->py())/cand->pt() * cand->pz()/cand->pt())>max_Dz_) continue;
139 
140  // dxy significance cut (safeguard against bizarre values)
141  if (min_DxySig_ > 0 && (tk->dxyError() <= 0 || std::abs(tk->dxy(beamSpot.position())/tk->dxyError()) < min_DxySig_)) continue;
142 
143  //dxy beamspot cut
144  float absDxy = std::abs(tk->dxy(beamSpot.position()));
145  if (absDxy > max_DXYBeamSpot_ || absDxy < min_DXYBeamSpot_ ) continue;
146  }
147 
148  //min muon hits cut
149  const reco::HitPattern& trackHits = tk->hitPattern();
150  if (trackHits.numberOfValidMuonHits() < min_NmuonHits_ ) continue;
151 
152  //pt difference cut
153  double candPt = cand->pt();
154  double trackPt = tk->pt();
155 
156  if (std::abs(candPt - trackPt) > max_PtDifference_ ) continue;
157 
158  //track pt cut
159  if (trackPt < min_TrackPt_ ) continue;
160 
161  // Pt threshold cut
162  double pt = cand->pt();
163  double err0 = tk->error(0);
164  double abspar0 = std::abs(tk->parameter(0));
165  double ptLx = pt;
166  // convert 50% efficiency threshold to 90% efficiency threshold
167  if (abspar0>0) ptLx += nsigma_Pt_*err0/abspar0*pt;
168  LogTrace("HLTMuonL3SimplePreFilter") << " ...Muon in loop, trackkRef pt= "
169  << tk->pt() << ", ptLx= " << ptLx
170  << " cand pT " << cand->pt();
171  if (ptLx<min_Pt_) continue;
172 
173  n++;
174  filterproduct.addObject(TriggerMuon,cand);
175  }//for iMu
176 
177  // filter decision
178  const bool accept (n >= min_N_);
179 
180  LogDebug("HLTMuonL3SimplePreFilter") << " >>>>> Result of HLTMuonL3PreFilter is " << accept << ", number of muons passing thresholds= " << n;
181 
182  return accept;
183 }
184 
185 
186 bool HLTMuonL3SimplePreFilter::triggerdByPreviousLevel(const reco::RecoChargedCandidateRef & candref, const std::vector<reco::RecoChargedCandidateRef>& vcands){
187  unsigned int i=0;
188  unsigned int i_max=vcands.size();
189  for (;i!=i_max;++i){
190  if (candref == vcands[i]) return true;
191  }
192 
193  return false;
194 }
195 
196 // declare this class as a framework plugin
#define LogDebug(id)
double z0() const
z coordinate
Definition: BeamSpot.h:68
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static bool triggerdByPreviousLevel(const reco::RecoChargedCandidateRef &, const std::vector< reco::RecoChargedCandidateRef > &)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
#define LogTrace(id)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > candToken_
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
HLTMuonL3SimplePreFilter(const edm::ParameterSet &)
~HLTMuonL3SimplePreFilter() override
double y0() const
y coordinate
Definition: BeamSpot.h:66
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const Point & position() const
position
Definition: BeamSpot.h:62
int numberOfValidMuonHits() const
Definition: HitPattern.h:906
double x0() const
x coordinate
Definition: BeamSpot.h:64