1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
|
/** \class HLTMuonTrimuonL3Filter
*
* See header file for documentation
*
* \author J. Alcaraz, P. Garcia
*
*/
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
#include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
#include "DataFormats/MuonReco/interface/MuonFwd.h"
#include "HLTMuonTrimuonL3Filter.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
#include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
#include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DataFormats/Math/interface/deltaR.h"
using namespace edm;
using namespace std;
using namespace reco;
using namespace trigger;
//
// constructors and destructor
//
HLTMuonTrimuonL3Filter::HLTMuonTrimuonL3Filter(const edm::ParameterSet& iConfig)
: HLTFilter(iConfig),
idealMagneticFieldRecordToken_(esConsumes()),
beamspotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
beamspotToken_(consumes<reco::BeamSpot>(beamspotTag_)),
candTag_(iConfig.getParameter<edm::InputTag>("CandTag")),
candToken_(consumes<reco::RecoChargedCandidateCollection>(candTag_)),
previousCandTag_(iConfig.getParameter<InputTag>("PreviousCandTag")),
previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
fast_Accept_(iConfig.getParameter<bool>("FastAccept")),
max_Eta_(iConfig.getParameter<double>("MaxEta")),
min_Nhits_(iConfig.getParameter<int>("MinNhits")),
max_Dr_(iConfig.getParameter<double>("MaxDr")),
max_Dz_(iConfig.getParameter<double>("MaxDz")),
chargeOpt_(iConfig.getParameter<int>("ChargeOpt")),
min_PtTriplet_(iConfig.getParameter<double>("MinPtTriplet")),
min_PtMax_(iConfig.getParameter<double>("MinPtMax")),
min_PtMin_(iConfig.getParameter<double>("MinPtMin")),
min_InvMass_(iConfig.getParameter<double>("MinInvMass")),
max_InvMass_(iConfig.getParameter<double>("MaxInvMass")),
min_Acop_(iConfig.getParameter<double>("MinAcop")),
max_Acop_(iConfig.getParameter<double>("MaxAcop")),
min_PtBalance_(iConfig.getParameter<double>("MinPtBalance")),
max_PtBalance_(iConfig.getParameter<double>("MaxPtBalance")),
nsigma_Pt_(iConfig.getParameter<double>("NSigmaPt")),
max_DCAMuMu_(iConfig.getParameter<double>("MaxDCAMuMu")),
max_YTriplet_(iConfig.getParameter<double>("MaxRapidityTriplet")),
theL3LinksLabel(iConfig.getParameter<InputTag>("InputLinks")),
linkToken_(consumes<reco::MuonTrackLinksCollection>(theL3LinksLabel)) {
LogDebug("HLTMuonTrimuonL3Filter")
<< " CandTag/MinN/MaxEta/MinNhits/MaxDr/MaxDz/MinPt1/MinPt2/MinInvMass/MaxInvMass/MinAcop/MaxAcop/MinPtBalance/"
"MaxPtBalance/NSigmaPt/MaxDzMuMu/MaxRapidityTriplet : "
<< candTag_.encode() << " " << fast_Accept_ << " " << max_Eta_ << " " << min_Nhits_ << " " << max_Dr_ << " "
<< max_Dz_ << " " << chargeOpt_ << " " << min_PtTriplet_ << " " << min_PtMax_ << " " << min_PtMin_ << " "
<< min_InvMass_ << " " << max_InvMass_ << " " << min_Acop_ << " " << max_Acop_ << " " << min_PtBalance_ << " "
<< max_PtBalance_ << " " << nsigma_Pt_ << " " << max_DCAMuMu_ << " " << max_YTriplet_;
}
HLTMuonTrimuonL3Filter::~HLTMuonTrimuonL3Filter() = default;
void HLTMuonTrimuonL3Filter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
makeHLTFilterDescription(desc);
desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltOfflineBeamSpot"));
desc.add<edm::InputTag>("CandTag", edm::InputTag("hltL3MuonCandidates"));
desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
desc.add<bool>("FastAccept", false);
desc.add<double>("MaxEta", 2.5);
desc.add<int>("MinNhits", 0);
desc.add<double>("MaxDr", 2.0);
desc.add<double>("MaxDz", 9999.0);
desc.add<int>("ChargeOpt", 0);
desc.add<double>("MinPtTriplet", 0.0);
desc.add<double>("MinPtMax", 3.0);
desc.add<double>("MinPtMin", 3.0);
desc.add<double>("MinInvMass", 2.8);
desc.add<double>("MaxInvMass", 3.4);
desc.add<double>("MinAcop", -1.0);
desc.add<double>("MaxAcop", 3.15);
desc.add<double>("MinPtBalance", -1.0);
desc.add<double>("MaxPtBalance", 999999.0);
desc.add<double>("NSigmaPt", 0.0);
desc.add<double>("MaxDCAMuMu", 99999.9);
desc.add<double>("MaxRapidityTriplet", 999999.0);
desc.add<edm::InputTag>("InputLinks", edm::InputTag(""));
descriptions.add("hltMuonTrimuonL3Filter", desc);
}
//
// member functions
//
// ------------ method called to produce the data ------------
bool HLTMuonTrimuonL3Filter::hltFilter(edm::Event& iEvent,
const edm::EventSetup& iSetup,
trigger::TriggerFilterObjectWithRefs& filterproduct) const {
double const MuMass = 0.106;
double const MuMass2 = MuMass * MuMass;
// All HLT filters must create and fill an HLT filter object,
// recording any reconstructed physics objects satisfying (or not)
// this HLT filter, and place it in the Event.
// get hold of trks
Handle<RecoChargedCandidateCollection> mucands;
if (saveTags())
filterproduct.addCollectionTag(candTag_);
iEvent.getByToken(candToken_, mucands);
// Test to see if we can use L3MuonTrajectorySeeds:
if (mucands->empty())
return false;
auto const& tk = (*mucands)[0].track();
bool useL3MTS = false;
if (tk->seedRef().isNonnull()) {
auto a = dynamic_cast<const L3MuonTrajectorySeed*>(tk->seedRef().get());
useL3MTS = a != nullptr;
}
// sort them by L2Track
std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;
// If we can use L3MuonTrajectory seeds run the older code:
if (useL3MTS) {
unsigned int maxI = mucands->size();
for (unsigned int i = 0; i != maxI; i++) {
const TrackRef& tk = (*mucands)[i].track();
edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef =
tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
TrackRef staTrack = l3seedRef->l2Track();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands, i));
}
}
// Using normal TrajectorySeeds:
else {
// Read Links collection:
edm::Handle<reco::MuonTrackLinksCollection> links;
iEvent.getByToken(linkToken_, links);
// Loop over RecoChargedCandidates:
for (unsigned int i(0); i < mucands->size(); ++i) {
RecoChargedCandidateRef cand(mucands, i);
for (auto const& link : *links) {
TrackRef tk = cand->track();
// Using the same method that was used to create the links between L3 and L2
// ToDo: there should be a better way than dR,dPt matching
const reco::Track& globalTrack = *link.globalTrack();
float dR2 = deltaR2(tk->eta(), tk->phi(), globalTrack.eta(), globalTrack.phi());
float dPt = std::abs(tk->pt() - globalTrack.pt()) / tk->pt();
const TrackRef staTrack = link.standAloneTrack();
if (dR2 < 0.02 * 0.02 and dPt < 0.001) {
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
}
} //MTL loop
} //RCC loop
} //end of using normal TrajectorySeeds
Handle<TriggerFilterObjectWithRefs> previousLevelCands;
iEvent.getByToken(previousCandToken_, previousLevelCands);
BeamSpot beamSpot;
Handle<BeamSpot> recoBeamSpotHandle;
iEvent.getByToken(beamspotToken_, recoBeamSpotHandle);
beamSpot = *recoBeamSpotHandle;
// Needed for DCA calculation
auto const& bFieldHandle = iSetup.getHandle(idealMagneticFieldRecordToken_);
// needed to compare to L2
vector<RecoChargedCandidateRef> vl2cands;
previousLevelCands->getObjects(TriggerMuon, vl2cands);
// look at all mucands, check cuts and add to filter object
int n = 0;
double e1, e2, e3;
Particle::LorentzVector p, p1, p2, p3;
auto L2toL3s_it1 = L2toL3s.begin();
auto L2toL3s_end = L2toL3s.end();
bool atLeastOneTriplet = false;
for (; L2toL3s_it1 != L2toL3s_end; ++L2toL3s_it1) {
if (!triggeredByLevel2(L2toL3s_it1->first, vl2cands))
continue;
//loop over the L3Tk reconstructed for this L2.
unsigned int iTk1 = 0;
unsigned int maxItk1 = L2toL3s_it1->second.size();
for (; iTk1 != maxItk1; iTk1++) {
bool thisL3Index1isDone = false;
RecoChargedCandidateRef& cand1 = L2toL3s_it1->second[iTk1];
TrackRef tk1 = cand1->get<TrackRef>();
// eta cut
LogDebug("HLTMuonTrimuonL3Filter") << " 1st muon in loop: q*pt= " << tk1->charge() * tk1->pt() << " ("
<< cand1->charge() * cand1->pt() << ") "
<< ", eta= " << tk1->eta() << " (" << cand1->eta() << ") "
<< ", hits= " << tk1->numberOfValidHits();
if (fabs(cand1->eta()) > max_Eta_)
continue;
// cut on number of hits
if (tk1->numberOfValidHits() < min_Nhits_)
continue;
//dr cut
// if (fabs(tk1->d0())>max_Dr_) continue;
if (fabs((-(cand1->vx() - beamSpot.x0()) * cand1->py() + (cand1->vy() - beamSpot.y0()) * cand1->px()) /
cand1->pt()) > max_Dr_)
continue;
//dz cut
if (fabs((cand1->vz() - beamSpot.z0()) -
((cand1->vx() - beamSpot.x0()) * cand1->px() + (cand1->vy() - beamSpot.y0()) * cand1->py()) /
cand1->pt() * cand1->pz() / cand1->pt()) > max_Dz_)
continue;
// Pt threshold cut
double pt1 = cand1->pt();
// double err1 = tk1->error(0);
// double abspar1 = fabs(tk1->parameter(0));
double ptLx1 = pt1;
// Don't convert to 90% efficiency threshold
LogDebug("HLTMuonTrimuonL3Filter") << " ... 1st muon in loop, pt1= " << pt1 << ", ptLx1= " << ptLx1;
auto L2toL3s_it2 = L2toL3s_it1;
L2toL3s_it2++;
for (; L2toL3s_it2 != L2toL3s_end; ++L2toL3s_it2) {
if (!triggeredByLevel2(L2toL3s_it2->first, vl2cands))
continue;
//loop over the L3Tk reconstructed for this L2.
unsigned int iTk2 = 0;
unsigned int maxItk2 = L2toL3s_it2->second.size();
for (; iTk2 != maxItk2; iTk2++) {
RecoChargedCandidateRef& cand2 = L2toL3s_it2->second[iTk2];
TrackRef tk2 = cand2->get<TrackRef>();
// eta cut
LogDebug("HLTMuonTrimuonL3Filter") << " 2nd muon in loop: q*pt= " << tk2->charge() * tk2->pt() << " ("
<< cand2->charge() * cand2->pt() << ") "
<< ", eta= " << tk2->eta() << " (" << cand2->eta() << ") "
<< ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
if (fabs(cand2->eta()) > max_Eta_)
continue;
// cut on number of hits
if (tk2->numberOfValidHits() < min_Nhits_)
continue;
//dr cut
// if (fabs(tk2->d0())>max_Dr_) continue;
if (fabs((-(cand2->vx() - beamSpot.x0()) * cand2->py() + (cand2->vy() - beamSpot.y0()) * cand2->px()) /
cand2->pt()) > max_Dr_)
continue;
//dz cut
if (fabs((cand2->vz() - beamSpot.z0()) -
((cand2->vx() - beamSpot.x0()) * cand2->px() + (cand2->vy() - beamSpot.y0()) * cand2->py()) /
cand2->pt() * cand2->pz() / cand2->pt()) > max_Dz_)
continue;
// Pt threshold cut
double pt2 = cand2->pt();
// double err2 = tk2->error(0);
// double abspar2 = fabs(tk2->parameter(0));
double ptLx2 = pt2;
// Don't convert to 90% efficiency threshold
LogDebug("HLTMuonTrimuonL3Filter") << " ... 2nd muon in loop, pt2= " << pt2 << ", ptLx2= " << ptLx2;
auto L2toL3s_it3 = L2toL3s_it2;
L2toL3s_it3++;
for (; L2toL3s_it3 != L2toL3s_end; ++L2toL3s_it3) {
if (!triggeredByLevel2(L2toL3s_it3->first, vl2cands))
continue;
//loop over the L3Tk reconstructed for this L2.
unsigned int iTk3 = 0;
unsigned int maxItk3 = L2toL3s_it3->second.size();
for (; iTk3 != maxItk3; iTk3++) {
RecoChargedCandidateRef& cand3 = L2toL3s_it3->second[iTk3];
TrackRef tk3 = cand3->get<TrackRef>();
// eta cut
LogDebug("HLTMuonTrimuonL3Filter") << " 3rd muon in loop: q*pt= " << tk3->charge() * tk3->pt() << " ("
<< cand3->charge() * cand3->pt() << ") "
<< ", eta= " << tk3->eta() << " (" << cand3->eta() << ") "
<< ", hits= " << tk3->numberOfValidHits();
if (fabs(cand3->eta()) > max_Eta_)
continue;
// cut on number of hits
if (tk3->numberOfValidHits() < min_Nhits_)
continue;
//dr cut
// if (fabs(tk1->d0())>max_Dr_) continue;
if (fabs((-(cand3->vx() - beamSpot.x0()) * cand3->py() + (cand3->vy() - beamSpot.y0()) * cand3->px()) /
cand3->pt()) > max_Dr_)
continue;
//dz cut
if (fabs((cand3->vz() - beamSpot.z0()) -
((cand3->vx() - beamSpot.x0()) * cand3->px() + (cand3->vy() - beamSpot.y0()) * cand3->py()) /
cand3->pt() * cand3->pz() / cand3->pt()) > max_Dz_)
continue;
// Pt threshold cut
double pt3 = cand3->pt();
// double err3 = tk3->error(0);
// double abspar3 = fabs(tk3->parameter(0));
double ptLx3 = pt3;
// Don't convert to 90% efficiency threshold
LogDebug("HLTMuonTrimuonL3Filter") << " ... 3rd muon in loop, pt3= " << pt3 << ", ptLx3= " << ptLx3;
if (ptLx1 > ptLx2 && ptLx1 > ptLx3 && ptLx1 < min_PtMax_)
continue;
else if (ptLx2 > ptLx1 && ptLx2 > ptLx3 && ptLx2 < min_PtMax_)
continue;
else if (ptLx3 < ptLx2 && ptLx3 > ptLx1 && ptLx3 < min_PtMax_)
continue;
if (ptLx1 < ptLx2 && ptLx1 < ptLx3 && ptLx1 < min_PtMin_)
continue;
else if (ptLx2 < ptLx1 && ptLx2 < ptLx3 && ptLx2 < min_PtMin_)
continue;
else if (ptLx3 < ptLx2 && ptLx3 < ptLx1 && ptLx3 < min_PtMin_)
continue;
if (chargeOpt_ > 0) {
if (abs(cand1->charge() + cand2->charge() + cand3->charge()) != chargeOpt_)
continue;
}
// Acoplanarity
double acop = fabs(cand1->phi() - cand2->phi());
if (acop > M_PI)
acop = 2 * M_PI - acop;
acop = M_PI - acop;
LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-2 acop= " << acop;
if (acop < min_Acop_)
continue;
if (acop > max_Acop_)
continue;
acop = fabs(cand1->phi() - cand3->phi());
if (acop > M_PI)
acop = 2 * M_PI - acop;
acop = M_PI - acop;
LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-3 acop= " << acop;
if (acop < min_Acop_)
continue;
if (acop > max_Acop_)
continue;
acop = fabs(cand3->phi() - cand2->phi());
if (acop > M_PI)
acop = 2 * M_PI - acop;
acop = M_PI - acop;
LogDebug("HLTMuonTrimuonL3Filter") << " ... 3-2 acop= " << acop;
if (acop < min_Acop_)
continue;
if (acop > max_Acop_)
continue;
// Pt balance
double ptbalance = fabs(cand1->pt() - cand2->pt());
if (ptbalance < min_PtBalance_)
continue;
if (ptbalance > max_PtBalance_)
continue;
ptbalance = fabs(cand1->pt() - cand3->pt());
if (ptbalance < min_PtBalance_)
continue;
if (ptbalance > max_PtBalance_)
continue;
ptbalance = fabs(cand3->pt() - cand2->pt());
if (ptbalance < min_PtBalance_)
continue;
if (ptbalance > max_PtBalance_)
continue;
// Combined trimuon system
e1 = sqrt(cand1->momentum().Mag2() + MuMass2);
e2 = sqrt(cand2->momentum().Mag2() + MuMass2);
e3 = sqrt(cand3->momentum().Mag2() + MuMass2);
p1 = Particle::LorentzVector(cand1->px(), cand1->py(), cand1->pz(), e1);
p2 = Particle::LorentzVector(cand2->px(), cand2->py(), cand2->pz(), e2);
p3 = Particle::LorentzVector(cand3->px(), cand3->py(), cand3->pz(), e3);
p = p1 + p2 + p3;
double pt123 = p.pt();
LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-2 pt123= " << pt123;
if (pt123 < min_PtTriplet_)
continue;
double invmass = abs(p.mass());
// if (invmass>0) invmass = sqrt(invmass); else invmass = 0;
LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-2 invmass= " << invmass;
if (invmass < min_InvMass_)
continue;
if (invmass > max_InvMass_)
continue;
// Delta Z between the two muons
//double DeltaZMuMu = fabs(tk2->dz(beamSpot.position())-tk1->dz(beamSpot.position()));
//if ( DeltaZMuMu > max_DzMuMu_) continue;
// DCA between the three muons
TransientTrack mu1TT(*tk1, &(*bFieldHandle));
TransientTrack mu2TT(*tk2, &(*bFieldHandle));
TransientTrack mu3TT(*tk3, &(*bFieldHandle));
TrajectoryStateClosestToPoint mu1TS = mu1TT.impactPointTSCP();
TrajectoryStateClosestToPoint mu2TS = mu2TT.impactPointTSCP();
TrajectoryStateClosestToPoint mu3TS = mu3TT.impactPointTSCP();
if (mu1TS.isValid() && mu2TS.isValid() && mu3TS.isValid()) {
ClosestApproachInRPhi cApp;
cApp.calculate(mu1TS.theState(), mu2TS.theState());
if (!cApp.status() || cApp.distance() > max_DCAMuMu_)
continue;
cApp.calculate(mu1TS.theState(), mu3TS.theState());
if (!cApp.status() || cApp.distance() > max_DCAMuMu_)
continue;
cApp.calculate(mu3TS.theState(), mu2TS.theState());
if (!cApp.status() || cApp.distance() > max_DCAMuMu_)
continue;
}
// Max dimuon |rapidity|
double rapidity = fabs(p.Rapidity());
if (rapidity > max_YTriplet_)
continue;
// Add this triplet
n++;
LogDebug("HLTMuonTrimuonL3Filter")
<< " Track1 passing filter: pt= " << cand1->pt() << ", eta: " << cand1->eta();
LogDebug("HLTMuonTrimuonL3Filter")
<< " Track2 passing filter: pt= " << cand2->pt() << ", eta: " << cand2->eta();
LogDebug("HLTMuonTrimuonL3Filter")
<< " Track2 passing filter: pt= " << cand3->pt() << ", eta: " << cand3->eta();
LogDebug("HLTMuonTrimuonL3Filter") << " Invmass= " << invmass;
bool i1done = false;
bool i2done = false;
bool i3done = false;
vector<RecoChargedCandidateRef> vref;
filterproduct.getObjects(TriggerMuon, vref);
for (auto& i : vref) {
RecoChargedCandidateRef candref = RecoChargedCandidateRef(i);
TrackRef tktmp = candref->get<TrackRef>();
if (tktmp == tk1) {
i1done = true;
} else if (tktmp == tk2) {
i2done = true;
} else if (tktmp == tk3) {
i3done = true;
}
if (i1done && i2done && i3done)
break;
}
if (!i1done) {
filterproduct.addObject(TriggerMuon, cand1);
}
if (!i2done) {
filterproduct.addObject(TriggerMuon, cand2);
}
if (!i3done) {
filterproduct.addObject(TriggerMuon, cand3);
}
//break anyway since a L3 track triplet has been found matching the criteria
thisL3Index1isDone = true;
atLeastOneTriplet = true;
break;
} //loop on the track of the third L2
//break the loop if fast accept.
if (atLeastOneTriplet && fast_Accept_)
break;
} //loop on the third L2
} //loop on the track of the second L2
//break the loop if fast accept.
if (atLeastOneTriplet && fast_Accept_)
break;
} //loop on the second L2
//break the loop if fast accept.
if (atLeastOneTriplet && fast_Accept_)
break;
if (thisL3Index1isDone)
break;
} //loop on tracks for first L2
//break the loop if fast accept.
if (atLeastOneTriplet && fast_Accept_)
break;
} //loop on the first L2
// filter decision
const bool accept(n >= 1);
LogDebug("HLTMuonTrimuonL3Filter") << " >>>>> Result of HLTMuonTrimuonL3Filter is " << accept
<< ", number of muon triplets passing thresholds= " << n;
return accept;
}
bool HLTMuonTrimuonL3Filter::triggeredByLevel2(const TrackRef& staTrack, vector<RecoChargedCandidateRef>& vcands) {
bool ok = false;
for (auto& vcand : vcands) {
if (vcand->get<TrackRef>() == staTrack) {
ok = true;
LogDebug("HLTMuonL3PreFilter") << "The L2 track triggered";
break;
}
}
return ok;
}
// declare this class as a framework plugin
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(HLTMuonTrimuonL3Filter);
|