CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTDisplacedmumumuVtxProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 
5 
13 
18 
24 
25 
27 
28 using namespace edm;
29 using namespace reco;
30 using namespace std;
31 using namespace trigger;
32 //
33 // constructors and destructor
34 //
36  src_ (iConfig.getParameter<edm::InputTag>("Src")),
37  previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
38  maxEta_ (iConfig.getParameter<double>("MaxEta")),
39  minPt_ (iConfig.getParameter<double>("MinPt")),
40  minPtTriplet_ (iConfig.getParameter<double>("MinPtTriplet")),
41  minInvMass_ (iConfig.getParameter<double>("MinInvMass")),
42  maxInvMass_ (iConfig.getParameter<double>("MaxInvMass")),
43  chargeOpt_ (iConfig.getParameter<int>("ChargeOpt"))
44 {
45  produces<VertexCollection>();
46 }
47 
48 
50 {
51 
52 }
53 
54 
55 // ------------ method called once each job just before starting event loop ------------
57 {
58 
59 }
60 
61 // ------------ method called once each job just after ending the event loop ------------
63 {
64 
65 }
66 
67 // ------------ method called on each new Event ------------
69 {
70  double const MuMass = 0.106;
71  double const MuMass2 = MuMass*MuMass;
72 
73 
74  // get hold of muon trks
76  iEvent.getByLabel (src_,mucands);
77 
78  //get the transient track builder:
80  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
81 
82  std::auto_ptr<VertexCollection> vertexCollection(new VertexCollection());
83 
84  // look at all mucands, check cuts and make vertices
85  double e1,e2,e3;
87 
88  RecoChargedCandidateCollection::const_iterator cand1;
89  RecoChargedCandidateCollection::const_iterator cand2;
90  RecoChargedCandidateCollection::const_iterator cand3;
91 
92  // get the objects passing the previous filter
94  iEvent.getByLabel (previousCandTag_,previousCands);
95 
96  vector<RecoChargedCandidateRef> vPrevCands;
97  previousCands->getObjects(TriggerMuon,vPrevCands);
98 
99  for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
100  TrackRef tk1 = cand1->get<TrackRef>();
101  LogDebug("HLTDisplacedMumumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
102 
103  //first check if this muon passed the previous filter
104  if( ! checkPreviousCand( tk1, vPrevCands) ) continue;
105 
106  // cuts
107  if (fabs(cand1->eta())>maxEta_) continue;
108  if (cand1->pt() < minPt_) continue;
109 
110  cand2 = cand1; cand2++;
111  for (; cand2!=mucands->end(); cand2++) {
112  TrackRef tk2 = cand2->get<TrackRef>();
113 
114  // eta cut
115  LogDebug("HLTMuonDimuonFilter") << " 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
116  //first check if this muon passed the previous filter
117  if( ! checkPreviousCand( tk2, vPrevCands) ) continue;
118 
119  // cuts
120  if (fabs(cand2->eta())>maxEta_) continue;
121  if (cand2->pt() < minPt_) continue;
122 
123  cand3 = cand2; cand3++;
124  for (; cand3!=mucands->end(); cand3++) {
125  TrackRef tk3 = cand3->get<TrackRef>();
126 
127  // eta cut
128  LogDebug("HLTMuonDimuonFilter") << " 3rd muon in loop: q*pt= " << cand3->charge()*cand3->pt() << ", eta= " << cand3->eta() << ", hits= " << tk3->numberOfValidHits() << ", d0= " << tk3->d0();
129  //first check if this muon passed the previous filter
130  if( ! checkPreviousCand( tk3, vPrevCands) ) continue;
131 
132  // cuts
133  if (fabs(cand3->eta())>maxEta_) continue;
134  if (cand3->pt() < minPt_) continue;
135 
136  // opposite sign or same sign
137  if (chargeOpt_>0) {
138  if (fabs (cand1->charge() + cand2->charge() + cand3->charge()) != chargeOpt_) continue;
139  }
140 
141  // Combined dimuon system
142  e1 = sqrt(cand1->momentum().Mag2()+MuMass2);
143  e2 = sqrt(cand2->momentum().Mag2()+MuMass2);
144  e3 = sqrt(cand3->momentum().Mag2()+MuMass2);
145  p1 = Particle::LorentzVector(cand1->px(),cand1->py(),cand1->pz(),e1);
146  p2 = Particle::LorentzVector(cand2->px(),cand2->py(),cand2->pz(),e2);
147  p3 = Particle::LorentzVector(cand3->px(),cand3->py(),cand3->pz(),e3);
148  p = p1+p2+p3;
149 
150 
151  if (p.pt()<minPtTriplet_) continue;
152 
153  double invmass = abs(p.mass());
154  LogDebug("HLTDisplacedMumumuFilter") << " ... 1-2 invmass= " << invmass;
155 
156  if (invmass<minInvMass_) continue;
157  if (invmass>maxInvMass_) continue;
158 
159  // do the vertex fit
160  vector<TransientTrack> t_tks;
161  TransientTrack ttkp1 = (*theB).build(&tk1);
162  TransientTrack ttkp2 = (*theB).build(&tk2);
163  TransientTrack ttkp3 = (*theB).build(&tk3);
164  t_tks.push_back(ttkp1);
165  t_tks.push_back(ttkp2);
166  t_tks.push_back(ttkp3);
167 
168 
169  if (t_tks.size()!=3) continue;
170 
171  KalmanVertexFitter kvf;
172  TransientVertex tv = kvf.vertex(t_tks);
173 
174  if (!tv.isValid()) continue;
175 
176  Vertex vertex = tv;
177 
178  // put vertex in the event
179  vertexCollection->push_back(vertex);
180  }
181  }
182  }
183  iEvent.put(vertexCollection);
184 }
185 
186 
187 
188 bool HLTDisplacedmumumuVtxProducer::checkPreviousCand(const TrackRef& trackref, vector<RecoChargedCandidateRef> & refVect){
189  bool ok=false;
190  for (unsigned int i=0; i<refVect.size(); i++) {
191  if ( refVect[i]->get<TrackRef>() == trackref ) {
192  ok=true;
193  break;
194  }
195  }
196  return ok;
197 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
HLTDisplacedmumumuVtxProducer(const edm::ParameterSet &)
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
#define abs(x)
Definition: mlp_lapack.h:159
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
tuple vertexCollection
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
T sqrt(T t)
Definition: SSEVec.h:46
bool checkPreviousCand(const reco::TrackRef &trackref, std::vector< reco::RecoChargedCandidateRef > &ref2)
double p2[4]
Definition: TauolaWrapper.h:90
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
const T & get() const
Definition: EventSetup.h:55
double p1[4]
Definition: TauolaWrapper.h:89
virtual void produce(edm::Event &, const edm::EventSetup &)
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:242
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
bool isValid() const
math::PtEtaPhiELorentzVectorF LorentzVector
double p3[4]
Definition: TauolaWrapper.h:91