CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTDisplacedmumuVtxProducer.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  minPtPair_ (iConfig.getParameter<double>("MinPtPair")),
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;
87 
88  RecoChargedCandidateCollection::const_iterator cand1;
89  RecoChargedCandidateCollection::const_iterator cand2;
90 
91  // get the objects passing the previous filter
93  iEvent.getByLabel (previousCandTag_,previousCands);
94 
95  vector<RecoChargedCandidateRef> vPrevCands;
96  previousCands->getObjects(TriggerMuon,vPrevCands);
97 
98  for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
99  TrackRef tk1 = cand1->get<TrackRef>();
100  LogDebug("HLTDisplacedMumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
101 
102  //first check if this muon passed the previous filter
103  if( ! checkPreviousCand( tk1, vPrevCands) ) continue;
104 
105  // cuts
106  if (fabs(cand1->eta())>maxEta_) continue;
107  if (cand1->pt() < minPt_) continue;
108 
109  cand2 = cand1; cand2++;
110  for (; cand2!=mucands->end(); cand2++) {
111  TrackRef tk2 = cand2->get<TrackRef>();
112 
113  // eta cut
114  LogDebug("HLTDisplacedmumuVtxProducer") << " 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
115  //first check if this muon passed the previous filter
116  if( ! checkPreviousCand( tk2, vPrevCands) ) continue;
117 
118  // cuts
119  if (fabs(cand2->eta())>maxEta_) continue;
120  if (cand2->pt() < minPt_) continue;
121 
122  // opposite sign or same sign
123  if (chargeOpt_<0) {
124  if (cand1->charge()*cand2->charge()>0) continue;
125  } else if (chargeOpt_>0) {
126  if (cand1->charge()*cand2->charge()<0) continue;
127  }
128 
129  // Combined dimuon system
130  e1 = sqrt(cand1->momentum().Mag2()+MuMass2);
131  e2 = sqrt(cand2->momentum().Mag2()+MuMass2);
132  p1 = Particle::LorentzVector(cand1->px(),cand1->py(),cand1->pz(),e1);
133  p2 = Particle::LorentzVector(cand2->px(),cand2->py(),cand2->pz(),e2);
134  p = p1+p2;
135 
136 
137  if (p.pt()<minPtPair_) continue;
138 
139  double invmass = abs(p.mass());
140  LogDebug("HLTDisplacedMumuFilter") << " ... 1-2 invmass= " << invmass;
141 
142  if (invmass<minInvMass_) continue;
143  if (invmass>maxInvMass_) continue;
144 
145  // do the vertex fit
146  vector<TransientTrack> t_tks;
147  TransientTrack ttkp1 = (*theB).build(&tk1);
148  TransientTrack ttkp2 = (*theB).build(&tk2);
149  t_tks.push_back(ttkp1);
150  t_tks.push_back(ttkp2);
151 
152 
153  if (t_tks.size()!=2) continue;
154 
155  KalmanVertexFitter kvf;
156  TransientVertex tv = kvf.vertex(t_tks);
157 
158  if (!tv.isValid()) continue;
159 
160  Vertex vertex = tv;
161 
162  // put vertex in the event
163  vertexCollection->push_back(vertex);
164  }
165  }
166  iEvent.put(vertexCollection);
167 }
168 
169 
170 
171 bool HLTDisplacedmumuVtxProducer::checkPreviousCand(const TrackRef& trackref, vector<RecoChargedCandidateRef> & refVect){
172  bool ok=false;
173  for (unsigned int i=0; i<refVect.size(); i++) {
174  if ( refVect[i]->get<TrackRef>() == trackref ) {
175  ok=true;
176  break;
177  }
178  }
179  return ok;
180 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
bool checkPreviousCand(const reco::TrackRef &trackref, std::vector< reco::RecoChargedCandidateRef > &ref2)
HLTDisplacedmumuVtxProducer(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
virtual void produce(edm::Event &, const edm::EventSetup &)
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
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