CMS 3D CMS Logo

HLTDisplacedmumuVtxProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 
6 
14 
19 
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  : srcTag_(iConfig.getParameter<edm::InputTag>("Src")),
37  srcToken_(consumes<reco::RecoChargedCandidateCollection>(srcTag_)),
38  previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
39  previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
40  maxEta_(iConfig.getParameter<double>("MaxEta")),
41  minPt_(iConfig.getParameter<double>("MinPt")),
42  minPtPair_(iConfig.getParameter<double>("MinPtPair")),
43  minInvMass_(iConfig.getParameter<double>("MinInvMass")),
44  maxInvMass_(iConfig.getParameter<double>("MaxInvMass")),
45  chargeOpt_(iConfig.getParameter<int>("ChargeOpt")) {
46  produces<VertexCollection>();
47 }
48 
50 
53  desc.add<edm::InputTag>("Src", edm::InputTag("hltL3MuonCandidates"));
54  desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
55  desc.add<double>("MaxEta", 2.5);
56  desc.add<double>("MinPt", 0.0);
57  desc.add<double>("MinPtPair", 0.0);
58  desc.add<double>("MinInvMass", 1.0);
59  desc.add<double>("MaxInvMass", 20.0);
60  desc.add<int>("ChargeOpt", -1);
61  descriptions.add("hltDisplacedmumuVtxProducer", desc);
62 }
63 
64 // ------------ method called on each new Event ------------
66  double const MuMass = 0.106;
67  double const MuMass2 = MuMass * MuMass;
68 
69  // get hold of muon trks
71  iEvent.getByToken(srcToken_, mucands);
72 
73  //get the transient track builder:
75  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
76 
77  std::unique_ptr<VertexCollection> vertexCollection(new VertexCollection());
78 
79  // look at all mucands, check cuts and make vertices
80  double e1, e2;
82 
83  RecoChargedCandidateCollection::const_iterator cand1;
84  RecoChargedCandidateCollection::const_iterator cand2;
85 
86  // get the objects passing the previous filter
88  iEvent.getByToken(previousCandToken_, previousCands);
89 
90  vector<RecoChargedCandidateRef> vPrevCands;
91  previousCands->getObjects(TriggerMuon, vPrevCands);
92 
93  for (cand1 = mucands->begin(); cand1 != mucands->end(); cand1++) {
94  TrackRef tk1 = cand1->get<TrackRef>();
95  LogDebug("HLTDisplacedMumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge() * cand1->pt()
96  << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
97 
98  //first check if this muon passed the previous filter
99  if (!checkPreviousCand(tk1, vPrevCands))
100  continue;
101 
102  // cuts
103  if (fabs(cand1->eta()) > maxEta_)
104  continue;
105  if (cand1->pt() < minPt_)
106  continue;
107 
108  cand2 = cand1;
109  cand2++;
110  for (; cand2 != mucands->end(); cand2++) {
111  TrackRef tk2 = cand2->get<TrackRef>();
112 
113  // eta cut
114  LogDebug("HLTDisplacedmumuVtxProducer")
115  << " 2nd muon in loop: q*pt= " << cand2->charge() * cand2->pt() << ", eta= " << cand2->eta()
116  << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
117  //first check if this muon passed the previous filter
118  if (!checkPreviousCand(tk2, vPrevCands))
119  continue;
120 
121  // cuts
122  if (fabs(cand2->eta()) > maxEta_)
123  continue;
124  if (cand2->pt() < minPt_)
125  continue;
126 
127  // opposite sign or same sign
128  if (chargeOpt_ < 0) {
129  if (cand1->charge() * cand2->charge() > 0)
130  continue;
131  } else if (chargeOpt_ > 0) {
132  if (cand1->charge() * cand2->charge() < 0)
133  continue;
134  }
135 
136  // Combined dimuon system
137  e1 = sqrt(cand1->momentum().Mag2() + MuMass2);
138  e2 = sqrt(cand2->momentum().Mag2() + MuMass2);
139  p1 = Particle::LorentzVector(cand1->px(), cand1->py(), cand1->pz(), e1);
140  p2 = Particle::LorentzVector(cand2->px(), cand2->py(), cand2->pz(), e2);
141  p = p1 + p2;
142 
143  if (p.pt() < minPtPair_)
144  continue;
145 
146  double invmass = abs(p.mass());
147  LogDebug("HLTDisplacedMumuFilter") << " ... 1-2 invmass= " << invmass;
148 
149  if (invmass < minInvMass_)
150  continue;
151  if (invmass > maxInvMass_)
152  continue;
153 
154  // do the vertex fit
155  vector<TransientTrack> t_tks;
156  TransientTrack ttkp1 = (*theB).build(&tk1);
157  TransientTrack ttkp2 = (*theB).build(&tk2);
158  t_tks.push_back(ttkp1);
159  t_tks.push_back(ttkp2);
160 
161  if (t_tks.size() != 2)
162  continue;
163 
164  KalmanVertexFitter kvf;
165  TransientVertex tv = kvf.vertex(t_tks);
166 
167  if (!tv.isValid())
168  continue;
169 
170  Vertex vertex = tv;
171 
172  // put vertex in the event
173  vertexCollection->push_back(vertex);
174  }
175  }
177 }
178 
180  const vector<RecoChargedCandidateRef>& refVect) const {
181  bool ok = false;
182  for (auto& i : refVect) {
183  if (i->get<TrackRef>() == trackref) {
184  ok = true;
185  break;
186  }
187  }
188  return ok;
189 }
ConfigurationDescriptions.h
edm::StreamID
Definition: StreamID.h:30
trigger::TriggerFilterObjectWithRefs
Definition: TriggerFilterObjectWithRefs.h:35
KalmanVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition: KalmanVertexFitter.h:49
mps_fire.i
i
Definition: mps_fire.py:355
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
MessageLogger.h
KalmanVertexFitter.h
ESHandle.h
TransientVertex::isValid
bool isValid() const
Definition: TransientVertex.h:195
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
HLTDisplacedmumuVtxProducer::minPtPair_
const double minPtPair_
Definition: HLTDisplacedmumuVtxProducer.h:46
HLTDisplacedmumuVtxProducer::maxInvMass_
const double maxInvMass_
Definition: HLTDisplacedmumuVtxProducer.h:48
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
edm::Ref::get
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
TriggerFilterObjectWithRefs.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
convertSQLiteXML.ok
bool ok
Definition: convertSQLiteXML.py:98
edm::Handle
Definition: AssociativeIterator.h:50
RecoCandidate.h
trigger::TriggerRefsCollections::getObjects
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
Definition: TriggerRefsCollections.h:452
reco::RecoChargedCandidateCollection
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
Definition: RecoChargedCandidateFwd.h:9
edm::Ref< TrackCollection >
HLTDisplacedmumuVtxProducer::srcToken_
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > srcToken_
Definition: HLTDisplacedmumuVtxProducer.h:41
CandidateFwd.h
HLTDisplacedmumuVtxProducer::chargeOpt_
const int chargeOpt_
Definition: HLTDisplacedmumuVtxProducer.h:49
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
TrackFwd.h
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
TransientTrackRecord
Definition: TransientTrackRecord.h:11
trigger::TriggerMuon
Definition: TriggerTypeDefs.h:68
edm::ESHandle< TransientTrackBuilder >
p2
double p2[4]
Definition: TauolaWrapper.h:90
HLTDisplacedmumuVtxProducer::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: HLTDisplacedmumuVtxProducer.cc:65
RefToBase.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
Vertex.h
bphysicsOniaDQM_cfi.vertex
vertex
Definition: bphysicsOniaDQM_cfi.py:7
TransientTrackBuilder.h
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
edm::ParameterSet
Definition: ParameterSet.h:36
HLTDisplacedmumuVtxProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTDisplacedmumuVtxProducer.cc:51
StorageManager_cfg.e1
e1
Definition: StorageManager_cfg.py:16
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
p1
double p1[4]
Definition: TauolaWrapper.h:89
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
TransientTrackRecord.h
get
#define get
HLTDisplacedmumuVtxProducer::HLTDisplacedmumuVtxProducer
HLTDisplacedmumuVtxProducer(const edm::ParameterSet &)
Definition: HLTDisplacedmumuVtxProducer.cc:35
HLTDisplacedmumuVtxProducer::maxEta_
const double maxEta_
Definition: HLTDisplacedmumuVtxProducer.h:44
reco::JetExtendedAssociation::LorentzVector
math::PtEtaPhiELorentzVectorF LorentzVector
Definition: JetExtendedAssociation.h:25
VertexFwd.h
TransientVertex.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
reco::TransientTrack
Definition: TransientTrack.h:19
RecoChargedCandidate.h
spclusmultinvestigator_cfi.vertexCollection
vertexCollection
Definition: spclusmultinvestigator_cfi.py:4
TriggerRefsCollections.h
HLTDisplacedmumuVtxProducer.h
trigger
Definition: HLTPrescaleTableCond.h:8
HLTDisplacedmumuVtxProducer::minInvMass_
const double minInvMass_
Definition: HLTDisplacedmumuVtxProducer.h:47
Candidate.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Event
Definition: Event.h:73
RecoChargedCandidateFwd.h
HLTDisplacedmumuVtxProducer::~HLTDisplacedmumuVtxProducer
~HLTDisplacedmumuVtxProducer() override
edm::InputTag
Definition: InputTag.h:15
reco::Vertex
Definition: Vertex.h:35
HLTDisplacedmumuVtxProducer::checkPreviousCand
bool checkPreviousCand(const reco::TrackRef &trackref, const std::vector< reco::RecoChargedCandidateRef > &ref2) const
Definition: HLTDisplacedmumuVtxProducer.cc:179
HLTDisplacedmumuVtxProducer::minPt_
const double minPt_
Definition: HLTDisplacedmumuVtxProducer.h:45
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22
HLTDisplacedmumuVtxProducer::previousCandToken_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
Definition: HLTDisplacedmumuVtxProducer.h:43