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 
6 
14 
19 
25 
26 
28 
29 using namespace edm;
30 using namespace reco;
31 using namespace std;
32 using namespace trigger;
33 //
34 // constructors and destructor
35 //
37  srcTag_ (iConfig.getParameter<edm::InputTag>("Src")),
38  srcToken_(consumes<reco::RecoChargedCandidateCollection>(srcTag_)),
39  previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
40  previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
41  maxEta_ (iConfig.getParameter<double>("MaxEta")),
42  minPt_ (iConfig.getParameter<double>("MinPt")),
43  minPtTriplet_ (iConfig.getParameter<double>("MinPtTriplet")),
44  minInvMass_ (iConfig.getParameter<double>("MinInvMass")),
45  maxInvMass_ (iConfig.getParameter<double>("MaxInvMass")),
46  chargeOpt_ (iConfig.getParameter<int>("ChargeOpt"))
47 {
48  produces<VertexCollection>();
49 }
50 
51 
53 {
54 
55 }
56 
57 void
60  desc.add<edm::InputTag>("Src",edm::InputTag("hltL3MuonCandidates"));
61  desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag(""));
62  desc.add<double>("MaxEta",2.5);
63  desc.add<double>("MinPt",0.0);
64  desc.add<double>("MinPtTriplet",0.0);
65  desc.add<double>("MinInvMass",1.0);
66  desc.add<double>("MaxInvMass",20.0);
67  desc.add<int>("ChargeOpt",-1);
68  descriptions.add("hltDisplacedmumumuVtxProducer", desc);
69 }
70 
71 
72 // ------------ method called once each job just before starting event loop ------------
74 {
75 
76 }
77 
78 // ------------ method called once each job just after ending the event loop ------------
80 {
81 
82 }
83 
84 // ------------ method called on each new Event ------------
86 {
87  double const MuMass = 0.106;
88  double const MuMass2 = MuMass*MuMass;
89 
90 
91  // get hold of muon trks
93  iEvent.getByToken(srcToken_,mucands);
94 
95  //get the transient track builder:
97  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
98 
99  std::auto_ptr<VertexCollection> vertexCollection(new VertexCollection());
100 
101  // look at all mucands, check cuts and make vertices
102  double e1,e2,e3;
104 
105  RecoChargedCandidateCollection::const_iterator cand1;
106  RecoChargedCandidateCollection::const_iterator cand2;
107  RecoChargedCandidateCollection::const_iterator cand3;
108 
109  // get the objects passing the previous filter
111  iEvent.getByToken(previousCandToken_,previousCands);
112 
113  vector<RecoChargedCandidateRef> vPrevCands;
114  previousCands->getObjects(TriggerMuon,vPrevCands);
115 
116  for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
117  TrackRef tk1 = cand1->get<TrackRef>();
118  LogDebug("HLTDisplacedMumumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
119 
120  //first check if this muon passed the previous filter
121  if( ! checkPreviousCand( tk1, vPrevCands) ) continue;
122 
123  // cuts
124  if (fabs(cand1->eta())>maxEta_) continue;
125  if (cand1->pt() < minPt_) continue;
126 
127  cand2 = cand1; cand2++;
128  for (; cand2!=mucands->end(); cand2++) {
129  TrackRef tk2 = cand2->get<TrackRef>();
130 
131  // eta cut
132  LogDebug("HLTMuonDimuonFilter") << " 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
133  //first check if this muon passed the previous filter
134  if( ! checkPreviousCand( tk2, vPrevCands) ) continue;
135 
136  // cuts
137  if (fabs(cand2->eta())>maxEta_) continue;
138  if (cand2->pt() < minPt_) continue;
139 
140  cand3 = cand2; cand3++;
141  for (; cand3!=mucands->end(); cand3++) {
142  TrackRef tk3 = cand3->get<TrackRef>();
143 
144  // eta cut
145  LogDebug("HLTMuonDimuonFilter") << " 3rd muon in loop: q*pt= " << cand3->charge()*cand3->pt() << ", eta= " << cand3->eta() << ", hits= " << tk3->numberOfValidHits() << ", d0= " << tk3->d0();
146  //first check if this muon passed the previous filter
147  if( ! checkPreviousCand( tk3, vPrevCands) ) continue;
148 
149  // cuts
150  if (fabs(cand3->eta())>maxEta_) continue;
151  if (cand3->pt() < minPt_) continue;
152 
153  // opposite sign or same sign
154  if (chargeOpt_>0) {
155  if (fabs (cand1->charge() + cand2->charge() + cand3->charge()) != chargeOpt_) continue;
156  }
157 
158  // Combined dimuon system
159  e1 = sqrt(cand1->momentum().Mag2()+MuMass2);
160  e2 = sqrt(cand2->momentum().Mag2()+MuMass2);
161  e3 = sqrt(cand3->momentum().Mag2()+MuMass2);
162  p1 = Particle::LorentzVector(cand1->px(),cand1->py(),cand1->pz(),e1);
163  p2 = Particle::LorentzVector(cand2->px(),cand2->py(),cand2->pz(),e2);
164  p3 = Particle::LorentzVector(cand3->px(),cand3->py(),cand3->pz(),e3);
165  p = p1+p2+p3;
166 
167 
168  if (p.pt()<minPtTriplet_) continue;
169 
170  double invmass = abs(p.mass());
171  LogDebug("HLTDisplacedMumumuFilter") << " ... 1-2 invmass= " << invmass;
172 
173  if (invmass<minInvMass_) continue;
174  if (invmass>maxInvMass_) continue;
175 
176  // do the vertex fit
177  vector<TransientTrack> t_tks;
178  TransientTrack ttkp1 = (*theB).build(&tk1);
179  TransientTrack ttkp2 = (*theB).build(&tk2);
180  TransientTrack ttkp3 = (*theB).build(&tk3);
181  t_tks.push_back(ttkp1);
182  t_tks.push_back(ttkp2);
183  t_tks.push_back(ttkp3);
184 
185 
186  if (t_tks.size()!=3) continue;
187 
188  KalmanVertexFitter kvf;
189  TransientVertex tv = kvf.vertex(t_tks);
190 
191  if (!tv.isValid()) continue;
192 
193  Vertex vertex = tv;
194 
195  // put vertex in the event
196  vertexCollection->push_back(vertex);
197  }
198  }
199  }
200  iEvent.put(vertexCollection);
201 }
202 
203 
204 
205 bool HLTDisplacedmumumuVtxProducer::checkPreviousCand(const TrackRef& trackref, vector<RecoChargedCandidateRef> & refVect){
206  bool ok=false;
207  for (unsigned int i=0; i<refVect.size(); i++) {
208  if ( refVect[i]->get<TrackRef>() == trackref ) {
209  ok=true;
210  break;
211  }
212  }
213  return ok;
214 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > srcToken_
HLTDisplacedmumumuVtxProducer(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
tuple vertexCollection
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool checkPreviousCand(const reco::TrackRef &trackref, std::vector< reco::RecoChargedCandidateRef > &ref2)
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:244
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double p2[4]
Definition: TauolaWrapper.h:90
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
const T & get() const
Definition: EventSetup.h:55
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double p1[4]
Definition: TauolaWrapper.h:89
virtual void produce(edm::Event &, const edm::EventSetup &)
bool isValid() const
math::PtEtaPhiELorentzVectorF LorentzVector
double p3[4]
Definition: TauolaWrapper.h:91