CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLTDisplacedmumuVtxProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 
21 
22 using namespace edm;
23 using namespace reco;
24 using namespace std;
25 using namespace trigger;
26 //
27 // constructors and destructor
28 //
30  : transientTrackRecordToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
31  srcTag_(iConfig.getParameter<edm::InputTag>("Src")),
32  srcToken_(consumes<reco::RecoChargedCandidateCollection>(srcTag_)),
33  previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
34  previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
35  matchToPrevious_(iConfig.getParameter<bool>("matchToPrevious")),
36  maxEta_(iConfig.getParameter<double>("MaxEta")),
37  minPt_(iConfig.getParameter<double>("MinPt")),
38  minPtPair_(iConfig.getParameter<double>("MinPtPair")),
39  minInvMass_(iConfig.getParameter<double>("MinInvMass")),
40  maxInvMass_(iConfig.getParameter<double>("MaxInvMass")),
41  chargeOpt_(iConfig.getParameter<int>("ChargeOpt")) {
42  produces<VertexCollection>();
43 }
44 
46 
49  desc.add<edm::InputTag>("Src", edm::InputTag("hltL3MuonCandidates"));
50  desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
51  desc.add<bool>("matchToPrevious", true);
52  desc.add<double>("MaxEta", 2.5);
53  desc.add<double>("MinPt", 0.0);
54  desc.add<double>("MinPtPair", 0.0);
55  desc.add<double>("MinInvMass", 1.0);
56  desc.add<double>("MaxInvMass", 20.0);
57  desc.add<int>("ChargeOpt", -1);
58  descriptions.add("hltDisplacedmumuVtxProducer", desc);
59 }
60 
61 // ------------ method called on each new Event ------------
63  double const MuMass = 0.106;
64  double const MuMass2 = MuMass * MuMass;
65 
66  // get hold of muon trks
68  iEvent.getByToken(srcToken_, mucands);
69 
70  // get the transient track builder
71  auto const& theB = iSetup.getHandle(transientTrackRecordToken_);
72 
73  std::unique_ptr<VertexCollection> vertexCollection(new VertexCollection());
74 
75  // look at all mucands, check cuts and make vertices
76  double e1, e2;
78 
79  RecoChargedCandidateCollection::const_iterator cand1;
80  RecoChargedCandidateCollection::const_iterator cand2;
81 
82  // get the objects passing the previous filter
84  if (matchToPrevious_)
85  iEvent.getByToken(previousCandToken_, previousCands);
86 
87  vector<RecoChargedCandidateRef> vPrevCands;
88  if (matchToPrevious_)
89  previousCands->getObjects(TriggerMuon, vPrevCands);
90 
91  for (cand1 = mucands->begin(); cand1 != mucands->end(); cand1++) {
92  TrackRef tk1 = cand1->get<TrackRef>();
93  LogDebug("HLTDisplacedMumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge() * cand1->pt()
94  << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
95 
96  //first check if this muon passed the previous filter
97  if (matchToPrevious_ && !checkPreviousCand(tk1, vPrevCands))
98  continue;
99 
100  // cuts
101  if (fabs(cand1->eta()) > maxEta_)
102  continue;
103  if (cand1->pt() < minPt_)
104  continue;
105 
106  cand2 = cand1;
107  cand2++;
108  for (; cand2 != mucands->end(); cand2++) {
109  TrackRef tk2 = cand2->get<TrackRef>();
110 
111  // eta cut
112  LogDebug("HLTDisplacedmumuVtxProducer")
113  << " 2nd muon in loop: q*pt= " << cand2->charge() * cand2->pt() << ", eta= " << cand2->eta()
114  << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
115  //first check if this muon passed the previous filter
116  if (matchToPrevious_ && !checkPreviousCand(tk2, vPrevCands))
117  continue;
118 
119  // cuts
120  if (fabs(cand2->eta()) > maxEta_)
121  continue;
122  if (cand2->pt() < minPt_)
123  continue;
124 
125  // opposite sign or same sign
126  if (chargeOpt_ < 0) {
127  if (cand1->charge() * cand2->charge() > 0)
128  continue;
129  } else if (chargeOpt_ > 0) {
130  if (cand1->charge() * cand2->charge() < 0)
131  continue;
132  }
133 
134  // Combined dimuon system
135  e1 = sqrt(cand1->momentum().Mag2() + MuMass2);
136  e2 = sqrt(cand2->momentum().Mag2() + MuMass2);
137  p1 = Particle::LorentzVector(cand1->px(), cand1->py(), cand1->pz(), e1);
138  p2 = Particle::LorentzVector(cand2->px(), cand2->py(), cand2->pz(), e2);
139  p = p1 + p2;
140 
141  if (p.pt() < minPtPair_)
142  continue;
143 
144  double invmass = abs(p.mass());
145  LogDebug("HLTDisplacedMumuFilter") << " ... 1-2 invmass= " << invmass;
146 
147  if (invmass < minInvMass_)
148  continue;
149  if (invmass > maxInvMass_)
150  continue;
151 
152  // do the vertex fit
153  vector<TransientTrack> t_tks;
154  TransientTrack ttkp1 = (*theB).build(&tk1);
155  TransientTrack ttkp2 = (*theB).build(&tk2);
156  t_tks.push_back(ttkp1);
157  t_tks.push_back(ttkp2);
158 
159  if (t_tks.size() != 2)
160  continue;
161 
162  KalmanVertexFitter kvf;
163  TransientVertex tv = kvf.vertex(t_tks);
164 
165  if (!tv.isValid())
166  continue;
167 
168  Vertex vertex = tv;
169 
170  // put vertex in the event
171  vertexCollection->push_back(vertex);
172  }
173  }
174  iEvent.put(std::move(vertexCollection));
175 }
176 
178  const vector<RecoChargedCandidateRef>& refVect) const {
179  bool ok = false;
180  for (auto& i : refVect) {
181  if (i->get<TrackRef>() == trackref) {
182  ok = true;
183  break;
184  }
185  }
186  return ok;
187 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const TString p2
Definition: fwPaths.cc:13
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
HLTDisplacedmumuVtxProducer(const edm::ParameterSet &)
std::vector< Vertex > VertexCollection
Definition: Vertex.h:12
tuple vertexCollection
bool checkPreviousCand(const reco::TrackRef &trackref, const std::vector< reco::RecoChargedCandidateRef > &ref2) const
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
def move
Definition: eostools.py:511
const TString p1
Definition: fwPaths.cc:12
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:232
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
~HLTDisplacedmumuVtxProducer() override
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackRecordToken_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
bool isValid() const
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
math::PtEtaPhiELorentzVectorF LorentzVector
#define LogDebug(id)