CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTmmkkFilter.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <cmath>
3 #include <vector>
4 
8 
11 
14 
23 
25 
27 
29 
30 using namespace edm;
31 using namespace reco;
32 using namespace std;
33 using namespace trigger;
34 
35 
36 // ----------------------------------------------------------------------
38  muCandTag_ (iConfig.getParameter<edm::InputTag>("MuCand")),
39  muCandToken_(consumes<reco::RecoChargedCandidateCollection>(muCandTag_)),
40  trkCandTag_ (iConfig.getParameter<edm::InputTag>("TrackCand")),
41  trkCandToken_(consumes<reco::RecoChargedCandidateCollection>(trkCandTag_)),
42  thirdTrackMass_(iConfig.getParameter<double>("ThirdTrackMass")),
43  fourthTrackMass_(iConfig.getParameter<double>("FourthTrackMass")),
44  maxEta_(iConfig.getParameter<double>("MaxEta")),
45  minPt_(iConfig.getParameter<double>("MinPt")),
46  minInvMass_(iConfig.getParameter<double>("MinInvMass")),
47  maxInvMass_(iConfig.getParameter<double>("MaxInvMass")),
48  maxNormalisedChi2_(iConfig.getParameter<double>("MaxNormalisedChi2")),
49  minLxySignificance_(iConfig.getParameter<double>("MinLxySignificance")),
50  minCosinePointingAngle_(iConfig.getParameter<double>("MinCosinePointingAngle")),
51  minD0Significance_(iConfig.getParameter<double>("MinD0Significance")),
52  fastAccept_(iConfig.getParameter<bool>("FastAccept")),
53  beamSpotTag_ (iConfig.getParameter<edm::InputTag> ("BeamSpotTag")),
54  beamSpotToken_(consumes<reco::BeamSpot>(beamSpotTag_))
55 {
56  produces<VertexCollection>();
57  produces<CandidateCollection>();
58 }
59 
60 
61 // ----------------------------------------------------------------------
63 
64 }
65 
66 void
70  desc.add<edm::InputTag>("MuCand",edm::InputTag("hltMuTracks"));
71  desc.add<edm::InputTag>("TrackCand",edm::InputTag("hltMumukAllConeTracks"));
72  desc.add<double>("ThirdTrackMass",0.106);
73  desc.add<double>("FourthTrackMass",0.106);
74  desc.add<double>("MaxEta",2.5);
75  desc.add<double>("MinPt",3.0);
76  desc.add<double>("MinInvMass",1.2);
77  desc.add<double>("MaxInvMass",2.2);
78  desc.add<double>("MaxNormalisedChi2",10.0);
79  desc.add<double>("MinLxySignificance",3.0);
80  desc.add<double>("MinCosinePointingAngle",0.9);
81  desc.add<double>("MinD0Significance",0.0);
82  desc.add<bool>("FastAccept",false);
83  desc.add<edm::InputTag>("BeamSpotTag",edm::InputTag("hltOfflineBeamSpot"));
84  descriptions.add("hltmmkkFilter",desc);
85 }
86 
87 // ----------------------------------------------------------------------
89 
90 }
91 
92 
93 // ----------------------------------------------------------------------
95 
96 }
97 
98 
99 // ----------------------------------------------------------------------
101 
102  const double MuMass(0.106);
103  const double MuMass2(MuMass*MuMass);
104 
105  const double thirdTrackMass2(thirdTrackMass_*thirdTrackMass_);
106  const double fourthTrackMass2(fourthTrackMass_*fourthTrackMass_);
107 
108  auto_ptr<CandidateCollection> output(new CandidateCollection());
109  auto_ptr<VertexCollection> vertexCollection(new VertexCollection());
110 
111  //get the transient track builder:
113  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
114 
115  //get the beamspot position
116  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
117  iEvent.getByToken(beamSpotToken_,recoBeamSpotHandle);
118  const reco::BeamSpot& vertexBeamSpot = *recoBeamSpotHandle;
119 
120  ESHandle<MagneticField> bFieldHandle;
121  iSetup.get<IdealMagneticFieldRecord>().get(bFieldHandle);
122 
123  const MagneticField* magField = bFieldHandle.product();
124 
125  TSCBLBuilderNoMaterial blsBuilder;
126 
127  // Ref to Candidate object to be recorded in filter object
130  RecoChargedCandidateRef refTrk;
131  RecoChargedCandidateRef refTrk2;
132 
133  // get hold of muon trks
135  iEvent.getByToken(muCandToken_,mucands);
136 
137  // get track candidates around displaced muons
139  iEvent.getByToken(trkCandToken_,trkcands);
140 
141  if (saveTags()) {
142  filterproduct.addCollectionTag(muCandTag_);
143  filterproduct.addCollectionTag(trkCandTag_);
144  }
145 
146  double e1,e2,e3,e4;
148 
149  //TrackRefs to mu cands in trkcand
150  vector<TrackRef> trkMuCands;
151 
152  //Already used mu tracks to avoid double counting candidates
153  vector<bool> isUsedCand(trkcands->size(),false);
154 
155  int counter = 0;
156 
157  //run algorithm
158  for (RecoChargedCandidateCollection::const_iterator mucand1=mucands->begin(), endCand1=mucands->end(); mucand1!=endCand1; ++mucand1) {
159 
160  if ( mucands->size()<2) break;
161  if ( trkcands->size()<2) break;
162 
163  TrackRef trk1 = mucand1->get<TrackRef>();
164  LogDebug("HLTDisplacedMumukkFilter") << " 1st muon: q*pt= " << trk1->charge()*trk1->pt() << ", eta= " << trk1->eta() << ", hits= " << trk1->numberOfValidHits();
165 
166  // eta cut
167  if (fabs(trk1->eta()) > maxEta_) continue;
168 
169  // Pt threshold cut
170  if (trk1->pt() < minPt_) continue;
171 
172  RecoChargedCandidateCollection::const_iterator mucand2 = mucand1; ++mucand2;
173 
174  for (RecoChargedCandidateCollection::const_iterator endCand2=mucands->end(); mucand2!=endCand2; ++mucand2) {
175 
176  TrackRef trk2 = mucand2->get<TrackRef>();
177 
178  LogDebug("HLTDisplacedMumukkFilter") << " 2nd muon: q*pt= " << trk2->charge()*trk2->pt() << ", eta= " << trk2->eta() << ", hits= " << trk2->numberOfValidHits();
179 
180  // eta cut
181  if (fabs(trk2->eta()) > maxEta_) continue;
182 
183  // Pt threshold cut
184  if (trk2->pt() < minPt_) continue;
185 
186  RecoChargedCandidateCollection::const_iterator trkcand, endCandTrk;
187 
188  std::vector<bool>::iterator isUsedIter, endIsUsedCand;
189 
190  //get overlapping muon candidates
191  for ( trkcand = trkcands->begin(), endCandTrk=trkcands->end(), isUsedIter = isUsedCand.begin(), endIsUsedCand = isUsedCand.end(); trkcand != endCandTrk && isUsedIter != endIsUsedCand; ++trkcand, ++isUsedIter) {
192  TrackRef trk3 = trkcand->get<TrackRef>();
193 
194  //check for overlapping muon tracks and save TrackRefs
195  if (overlap(*mucand1,*trkcand)) {
196  trkMuCands.push_back(trk3);
197  *isUsedIter = true;
198  continue;
199  }
200  else if (overlap(*mucand2,*trkcand)){
201  trkMuCands.push_back(trk3);
202  *isUsedIter = true;
203  continue;
204  }
205 
206  if(trkMuCands.size()==2) break;
207  }
208 
209  //Not all overlapping candidates found, skip event
210  //if (trkMuCands.size()!=2) continue;
211 
212  //combine muons with all tracks
213  for ( trkcand = trkcands->begin(), endCandTrk=trkcands->end(), isUsedIter = isUsedCand.begin(), endIsUsedCand = isUsedCand.end(); trkcand != endCandTrk && isUsedIter != endIsUsedCand; ++trkcand, ++isUsedIter) {
214 
215  TrackRef trk3 = trkcand->get<TrackRef>();
216 
217  LogDebug("HLTDisplacedMumukkFilter") << " 3rd track: q*pt= " << trk3->charge()*trk3->pt() << ", eta= " << trk3->eta() << ", hits= " << trk3->numberOfValidHits();
218 
219  //skip overlapping muon candidates
220  bool skip=false;
221  for (unsigned int itmc=0;itmc<trkMuCands.size();itmc++) if(trk3==trkMuCands.at(itmc)) skip=true;
222  if(skip) continue;
223 
224  //skip already used tracks
225  if(*isUsedIter) continue;
226 
227  // eta cut
228  if (fabs(trk3->eta()) > maxEta_) continue;
229 
230  // Pt threshold cut
231  if (trk3->pt() < minPt_) continue;
232 
233  RecoChargedCandidateCollection::const_iterator trkcand2;
234 
235  std::vector<bool>::iterator isUsedIter2;
236 
237  for ( trkcand2 = trkcands->begin(), isUsedIter2 = isUsedCand.begin() ; trkcand2 != endCandTrk && isUsedIter2 != endIsUsedCand; ++trkcand2,++isUsedIter2) {
238 
239  if (trkcand2==trkcand) continue;
240 
241  TrackRef trk4 = trkcand2->get<TrackRef>();
242 
243  LogDebug("HLTDisplacedMumukkFilter") << " 4th track: q*pt= " << trk4->charge()*trk4->pt() << ", eta= " << trk4->eta() << ", hits= " << trk4->numberOfValidHits();
244 
245  //skip overlapping muon candidates
246  bool skip2=false;
247  for (unsigned int itmc=0;itmc<trkMuCands.size();itmc++) if(trk4==trkMuCands.at(itmc)) skip2=true;
248  if(skip2) continue;
249 
250  //skip already used tracks
251  if(*isUsedIter2) continue;
252 
253  // eta cut
254  if (fabs(trk4->eta()) > maxEta_) continue;
255 
256  // Pt threshold cut
257  if (trk4->pt() < minPt_) continue;
258 
259  // Combined system
260  e1 = sqrt(trk1->momentum().Mag2()+MuMass2);
261  e2 = sqrt(trk2->momentum().Mag2()+MuMass2);
262  e3 = sqrt(trk3->momentum().Mag2()+thirdTrackMass2);
263  e4 = sqrt(trk4->momentum().Mag2()+fourthTrackMass2);
264 
265  p1 = Particle::LorentzVector(trk1->px(),trk1->py(),trk1->pz(),e1);
266  p2 = Particle::LorentzVector(trk2->px(),trk2->py(),trk2->pz(),e2);
267  p3 = Particle::LorentzVector(trk3->px(),trk3->py(),trk3->pz(),e3);
268  p4 = Particle::LorentzVector(trk4->px(),trk4->py(),trk4->pz(),e4);
269 
270  p = p1+p2+p3+p4;
271 
272  //invariant mass cut
273  double invmass = abs(p.mass());
274 
275  LogDebug("HLTDisplacedMumukkFilter") << " Invmass= " << invmass;
276 
277  FreeTrajectoryState InitialFTS = initialFreeState(*trk3, magField);
278  TrajectoryStateClosestToBeamLine tscb( blsBuilder(InitialFTS, *recoBeamSpotHandle) );
279  double d0sig = tscb.transverseImpactParameter().significance();
280 
281  if (d0sig < minD0Significance_) continue;
282 
283  FreeTrajectoryState InitialFTS2 = initialFreeState(*trk4, magField);
284  TrajectoryStateClosestToBeamLine tscb2( blsBuilder(InitialFTS2, *recoBeamSpotHandle) );
286 
287  if (d0sig < minD0Significance_) continue;
288 
289  if (invmass>maxInvMass_ || invmass<minInvMass_) continue;
290 
291  // do the vertex fit
292  vector<TransientTrack> t_tks;
293  t_tks.push_back((*theB).build(&trk1));
294  t_tks.push_back((*theB).build(&trk2));
295  t_tks.push_back((*theB).build(&trk3));
296  t_tks.push_back((*theB).build(&trk4));
297 
298  if (t_tks.size()!=4) continue;
299 
300  KalmanVertexFitter kvf;
301  TransientVertex tv = kvf.vertex(t_tks);
302 
303  if (!tv.isValid()) continue;
304 
305  Vertex vertex = tv;
306 
307  // get vertex position and error to calculate the decay length significance
308  GlobalPoint secondaryVertex = tv.position();
309  GlobalError err = tv.positionError();
310 
311  //calculate decay length significance w.r.t. the beamspot
312  GlobalPoint displacementFromBeamspot( -1*((vertexBeamSpot.x0() -secondaryVertex.x()) + (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dxdz()), -1*((vertexBeamSpot.y0() - secondaryVertex.y())+ (secondaryVertex.z() -vertexBeamSpot.z0()) * vertexBeamSpot.dydz()), 0);
313 
314  float lxy = displacementFromBeamspot.perp();
315  float lxyerr = sqrt(err.rerr(displacementFromBeamspot));
316 
317  // get normalizes chi2
318  float normChi2 = tv.normalisedChiSquared();
319 
320  //calculate the angle between the decay length and the mumu momentum
321  Vertex::Point vperp(displacementFromBeamspot.x(),displacementFromBeamspot.y(),0.);
322  math::XYZVector pperp(p.x(),p.y(),0.);
323 
324  float cosAlpha = vperp.Dot(pperp)/(vperp.R()*pperp.R());
325 
326  LogDebug("HLTDisplacedMumukkFilter") << " vertex fit normalised chi2: " << normChi2 << ", Lxy significance: " << lxy/lxyerr << ", cosine pointing angle: " << cosAlpha;
327 
328  // put vertex in the event
329  vertexCollection->push_back(vertex);
330 
331  if (normChi2 > maxNormalisedChi2_) continue;
332  if (lxy/lxyerr < minLxySignificance_) continue;
333  if(cosAlpha < minCosinePointingAngle_) continue;
334 
335  LogDebug("HLTDisplacedMumukkFilter") << " Event passed!";
336 
337  //Add event
338  ++counter;
339 
340  //Get refs
341  bool i1done = false;
342  bool i2done = false;
343  bool i3done = false;
344  bool i4done = false;
345  vector<RecoChargedCandidateRef> vref;
346  filterproduct.getObjects(TriggerMuon,vref);
347  for (unsigned int i=0; i<vref.size(); i++) {
349  TrackRef trktmp = candref->get<TrackRef>();
350  if (trktmp==trk1) {
351  i1done = true;
352  } else if (trktmp==trk2) {
353  i2done = true;
354  } else if (trktmp==trk3) {
355  i3done = true;
356  } else if (trktmp==trk4) {
357  i4done = true;
358  }
359  if (i1done && i2done && i3done && i4done) break;
360  }
361 
362  if (!i1done) {
363  refMu1=RecoChargedCandidateRef( Ref<RecoChargedCandidateCollection> (mucands,distance(mucands->begin(), mucand1)));
364  filterproduct.addObject(TriggerMuon,refMu1);
365  }
366  if (!i2done) {
367  refMu2=RecoChargedCandidateRef( Ref<RecoChargedCandidateCollection> (mucands,distance(mucands->begin(),mucand2)));
368  filterproduct.addObject(TriggerMuon,refMu2);
369  }
370  if (!i3done) {
371  refTrk=RecoChargedCandidateRef( Ref<RecoChargedCandidateCollection> (trkcands,distance(trkcands->begin(),trkcand)));
372  filterproduct.addObject(TriggerTrack,refTrk);
373  }
374  if (!i4done) {
375  refTrk2=RecoChargedCandidateRef( Ref<RecoChargedCandidateCollection> (trkcands,distance(trkcands->begin(),trkcand)));
376  filterproduct.addObject(TriggerTrack,refTrk2);
377  }
378 
379  if (fastAccept_) break;
380  }
381  }
382 
383  trkMuCands.clear();
384  }
385  }
386 
387  // filter decision
388  const bool accept (counter >= 1);
389 
390  LogDebug("HLTDisplacedMumukkFilter") << " >>>>> Result of HLTDisplacedMumukkFilter is "<< accept << ", number of muon pairs passing thresholds= " << counter;
391 
392  iEvent.put(vertexCollection);
393 
394  return accept;
395 }
396 
397 // ----------------------------------------------------------------------
399  const MagneticField* field)
400 {
401  Basic3DVector<float> pos( tk.vertex());
402  GlobalPoint gpos( pos);
403  Basic3DVector<float> mom( tk.momentum());
404  GlobalVector gmom( mom);
405  GlobalTrajectoryParameters par( gpos, gmom, tk.charge(), field);
407  return FreeTrajectoryState( par, err);
408 }
409 
411 
412  double eps(1.44e-4);
413 
414  double dpt = a.pt() - b.pt();
415  dpt *= dpt;
416 
417  double dphi = deltaPhi(a.phi(), b.phi());
418  dphi *= dphi;
419 
420  double deta = a.eta() - b.eta();
421  deta *= deta;
422 
423  if ((dphi + deta) < eps) {
424  return 1;
425  }
426 
427  return 0;
428 
429 }
430 
#define LogDebug(id)
GlobalError positionError() const
double z0() const
z coordinate
Definition: BeamSpot.h:68
int i
Definition: DBlmapReader.cc:9
const double maxInvMass_
Definition: HLTmmkkFilter.h:65
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
edm::InputTag muCandTag_
Definition: HLTmmkkFilter.h:55
T perp() const
Definition: PV3DBase.h:72
const double minLxySignificance_
Definition: HLTmmkkFilter.h:67
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual double pt() const =0
transverse momentum
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
T y() const
Definition: PV3DBase.h:63
edm::Ref< RecoChargedCandidateCollection > RecoChargedCandidateRef
reference to an object in a collection of RecoChargedCandidate objects
const double maxNormalisedChi2_
Definition: HLTmmkkFilter.h:66
virtual void endJob()
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:662
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
tuple vertexCollection
const double thirdTrackMass_
Definition: HLTmmkkFilter.h:60
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
static int overlap(const reco::Candidate &, const reco::Candidate &)
double dydz() const
dydz slope
Definition: BeamSpot.h:84
float normalisedChiSquared() const
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:674
int iEvent
Definition: GenABIO.cc:230
virtual void beginJob()
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > muCandToken_
Definition: HLTmmkkFilter.h:56
const bool fastAccept_
Definition: HLTmmkkFilter.h:70
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:718
T sqrt(T t)
Definition: SSEVec.h:48
double p4[4]
Definition: TauolaWrapper.h:92
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:64
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
Definition: HLTmmkkFilter.h:72
const double maxEta_
Definition: HLTmmkkFilter.h:62
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:244
edm::InputTag trkCandTag_
Definition: HLTmmkkFilter.h:57
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const double minInvMass_
Definition: HLTmmkkFilter.h:64
double p2[4]
Definition: TauolaWrapper.h:90
const double minPt_
Definition: HLTmmkkFilter.h:63
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
double significance() const
Definition: Measurement1D.h:32
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
const double minD0Significance_
Definition: HLTmmkkFilter.h:69
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
T rerr(const GlobalPoint &aPoint) const
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
double b
Definition: hdecay.h:120
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static FreeTrajectoryState initialFreeState(const reco::Track &, const MagneticField *)
const double minCosinePointingAngle_
Definition: HLTmmkkFilter.h:68
double value() const
Definition: Measurement1D.h:28
const double fourthTrackMass_
Definition: HLTmmkkFilter.h:61
bool saveTags() const
Definition: HLTFilter.h:45
double p1[4]
Definition: TauolaWrapper.h:89
HLTmmkkFilter(const edm::ParameterSet &)
double a
Definition: hdecay.h:121
static std::atomic< unsigned int > counter
double y0() const
y coordinate
Definition: BeamSpot.h:66
int charge() const
track electric charge
Definition: TrackBase.h:554
edm::OwnVector< Candidate > CandidateCollection
collection of Candidate objects
Definition: CandidateFwd.h:21
T x() const
Definition: PV3DBase.h:62
bool isValid() const
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > trkCandToken_
Definition: HLTmmkkFilter.h:58
math::PtEtaPhiELorentzVectorF LorentzVector
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
double p3[4]
Definition: TauolaWrapper.h:91
double x0() const
x coordinate
Definition: BeamSpot.h:64