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 
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  minPtPair_ (iConfig.getParameter<double>("MinPtPair")),
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>("MinPtPair",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("hltDisplacedmumuVtxProducer", desc);
69 }
70 
71 // ------------ method called once each job just before starting event loop ------------
73 {
74 
75 }
76 
77 // ------------ method called once each job just after ending the event loop ------------
79 {
80 
81 }
82 
83 // ------------ method called on each new Event ------------
85 {
86  double const MuMass = 0.106;
87  double const MuMass2 = MuMass*MuMass;
88 
89 
90  // get hold of muon trks
92  iEvent.getByToken(srcToken_,mucands);
93 
94  //get the transient track builder:
96  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
97 
98  std::auto_ptr<VertexCollection> vertexCollection(new VertexCollection());
99 
100  // look at all mucands, check cuts and make vertices
101  double e1,e2;
103 
104  RecoChargedCandidateCollection::const_iterator cand1;
105  RecoChargedCandidateCollection::const_iterator cand2;
106 
107  // get the objects passing the previous filter
109  iEvent.getByToken(previousCandToken_,previousCands);
110 
111  vector<RecoChargedCandidateRef> vPrevCands;
112  previousCands->getObjects(TriggerMuon,vPrevCands);
113 
114  for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
115  TrackRef tk1 = cand1->get<TrackRef>();
116  LogDebug("HLTDisplacedMumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
117 
118  //first check if this muon passed the previous filter
119  if( ! checkPreviousCand( tk1, vPrevCands) ) continue;
120 
121  // cuts
122  if (fabs(cand1->eta())>maxEta_) continue;
123  if (cand1->pt() < minPt_) continue;
124 
125  cand2 = cand1; cand2++;
126  for (; cand2!=mucands->end(); cand2++) {
127  TrackRef tk2 = cand2->get<TrackRef>();
128 
129  // eta cut
130  LogDebug("HLTDisplacedmumuVtxProducer") << " 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
131  //first check if this muon passed the previous filter
132  if( ! checkPreviousCand( tk2, vPrevCands) ) continue;
133 
134  // cuts
135  if (fabs(cand2->eta())>maxEta_) continue;
136  if (cand2->pt() < minPt_) continue;
137 
138  // opposite sign or same sign
139  if (chargeOpt_<0) {
140  if (cand1->charge()*cand2->charge()>0) continue;
141  } else if (chargeOpt_>0) {
142  if (cand1->charge()*cand2->charge()<0) continue;
143  }
144 
145  // Combined dimuon system
146  e1 = sqrt(cand1->momentum().Mag2()+MuMass2);
147  e2 = sqrt(cand2->momentum().Mag2()+MuMass2);
148  p1 = Particle::LorentzVector(cand1->px(),cand1->py(),cand1->pz(),e1);
149  p2 = Particle::LorentzVector(cand2->px(),cand2->py(),cand2->pz(),e2);
150  p = p1+p2;
151 
152 
153  if (p.pt()<minPtPair_) continue;
154 
155  double invmass = abs(p.mass());
156  LogDebug("HLTDisplacedMumuFilter") << " ... 1-2 invmass= " << invmass;
157 
158  if (invmass<minInvMass_) continue;
159  if (invmass>maxInvMass_) continue;
160 
161  // do the vertex fit
162  vector<TransientTrack> t_tks;
163  TransientTrack ttkp1 = (*theB).build(&tk1);
164  TransientTrack ttkp2 = (*theB).build(&tk2);
165  t_tks.push_back(ttkp1);
166  t_tks.push_back(ttkp2);
167 
168 
169  if (t_tks.size()!=2) 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  iEvent.put(vertexCollection);
183 }
184 
185 
186 
187 bool HLTDisplacedmumuVtxProducer::checkPreviousCand(const TrackRef& trackref, vector<RecoChargedCandidateRef> & refVect){
188  bool ok=false;
189  for (unsigned int i=0; i<refVect.size(); i++) {
190  if ( refVect[i]->get<TrackRef>() == trackref ) {
191  ok=true;
192  break;
193  }
194  }
195  return ok;
196 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
bool checkPreviousCand(const reco::TrackRef &trackref, std::vector< reco::RecoChargedCandidateRef > &ref2)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
HLTDisplacedmumuVtxProducer(const edm::ParameterSet &)
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
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
T sqrt(T t)
Definition: SSEVec.h:48
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > srcToken_
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:244
virtual void produce(edm::Event &, const edm::EventSetup &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double p2[4]
Definition: TauolaWrapper.h:90
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const T & get() const
Definition: EventSetup.h:56
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double p1[4]
Definition: TauolaWrapper.h:89
bool isValid() const
math::PtEtaPhiELorentzVectorF LorentzVector