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  matchToPrevious_(iConfig.getParameter<bool>("matchToPrevious")),
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  produces<VertexCollection>();
48 }
49 
51 
54  desc.add<edm::InputTag>("Src", edm::InputTag("hltL3MuonCandidates"));
55  desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
56  desc.add<bool>("matchToPrevious", true);
57  desc.add<double>("MaxEta", 2.5);
58  desc.add<double>("MinPt", 0.0);
59  desc.add<double>("MinPtPair", 0.0);
60  desc.add<double>("MinInvMass", 1.0);
61  desc.add<double>("MaxInvMass", 20.0);
62  desc.add<int>("ChargeOpt", -1);
63  descriptions.add("hltDisplacedmumuVtxProducer", desc);
64 }
65 
66 // ------------ method called on each new Event ------------
68  double const MuMass = 0.106;
69  double const MuMass2 = MuMass * MuMass;
70 
71  // get hold of muon trks
73  iEvent.getByToken(srcToken_, mucands);
74 
75  //get the transient track builder:
77  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
78 
79  std::unique_ptr<VertexCollection> vertexCollection(new VertexCollection());
80 
81  // look at all mucands, check cuts and make vertices
82  double e1, e2;
84 
85  RecoChargedCandidateCollection::const_iterator cand1;
86  RecoChargedCandidateCollection::const_iterator cand2;
87 
88  // get the objects passing the previous filter
90  if (matchToPrevious_)
91  iEvent.getByToken(previousCandToken_, previousCands);
92 
93  vector<RecoChargedCandidateRef> vPrevCands;
94  if (matchToPrevious_)
95  previousCands->getObjects(TriggerMuon, vPrevCands);
96 
97  for (cand1 = mucands->begin(); cand1 != mucands->end(); cand1++) {
98  TrackRef tk1 = cand1->get<TrackRef>();
99  LogDebug("HLTDisplacedMumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge() * cand1->pt()
100  << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
101 
102  //first check if this muon passed the previous filter
103  if (matchToPrevious_ && !checkPreviousCand(tk1, vPrevCands))
104  continue;
105 
106  // cuts
107  if (fabs(cand1->eta()) > maxEta_)
108  continue;
109  if (cand1->pt() < minPt_)
110  continue;
111 
112  cand2 = cand1;
113  cand2++;
114  for (; cand2 != mucands->end(); cand2++) {
115  TrackRef tk2 = cand2->get<TrackRef>();
116 
117  // eta cut
118  LogDebug("HLTDisplacedmumuVtxProducer")
119  << " 2nd muon in loop: q*pt= " << cand2->charge() * cand2->pt() << ", eta= " << cand2->eta()
120  << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
121  //first check if this muon passed the previous filter
122  if (matchToPrevious_ && !checkPreviousCand(tk2, vPrevCands))
123  continue;
124 
125  // cuts
126  if (fabs(cand2->eta()) > maxEta_)
127  continue;
128  if (cand2->pt() < minPt_)
129  continue;
130 
131  // opposite sign or same sign
132  if (chargeOpt_ < 0) {
133  if (cand1->charge() * cand2->charge() > 0)
134  continue;
135  } else if (chargeOpt_ > 0) {
136  if (cand1->charge() * cand2->charge() < 0)
137  continue;
138  }
139 
140  // Combined dimuon system
141  e1 = sqrt(cand1->momentum().Mag2() + MuMass2);
142  e2 = sqrt(cand2->momentum().Mag2() + MuMass2);
143  p1 = Particle::LorentzVector(cand1->px(), cand1->py(), cand1->pz(), e1);
144  p2 = Particle::LorentzVector(cand2->px(), cand2->py(), cand2->pz(), e2);
145  p = p1 + p2;
146 
147  if (p.pt() < minPtPair_)
148  continue;
149 
150  double invmass = abs(p.mass());
151  LogDebug("HLTDisplacedMumuFilter") << " ... 1-2 invmass= " << invmass;
152 
153  if (invmass < minInvMass_)
154  continue;
155  if (invmass > maxInvMass_)
156  continue;
157 
158  // do the vertex fit
159  vector<TransientTrack> t_tks;
160  TransientTrack ttkp1 = (*theB).build(&tk1);
161  TransientTrack ttkp2 = (*theB).build(&tk2);
162  t_tks.push_back(ttkp1);
163  t_tks.push_back(ttkp2);
164 
165  if (t_tks.size() != 2)
166  continue;
167 
168  KalmanVertexFitter kvf;
169  TransientVertex tv = kvf.vertex(t_tks);
170 
171  if (!tv.isValid())
172  continue;
173 
174  Vertex vertex = tv;
175 
176  // put vertex in the event
177  vertexCollection->push_back(vertex);
178  }
179  }
181 }
182 
184  const vector<RecoChargedCandidateRef>& refVect) const {
185  bool ok = false;
186  for (auto& i : refVect) {
187  if (i->get<TrackRef>() == trackref) {
188  ok = true;
189  break;
190  }
191  }
192  return ok;
193 }
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
electrons_cff.bool
bool
Definition: electrons_cff.py:393
mps_fire.i
i
Definition: mps_fire.py:428
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:47
HLTDisplacedmumuVtxProducer::maxInvMass_
const double maxInvMass_
Definition: HLTDisplacedmumuVtxProducer.h:49
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
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:50
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
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:67
RefToBase.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
Vertex.h
bphysicsOniaDQM_cfi.vertex
vertex
Definition: bphysicsOniaDQM_cfi.py:7
TransientTrackBuilder.h
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
edm::ParameterSet
Definition: ParameterSet.h:47
HLTDisplacedmumuVtxProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTDisplacedmumuVtxProducer.cc:52
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:45
reco::JetExtendedAssociation::LorentzVector
math::PtEtaPhiELorentzVectorF LorentzVector
Definition: JetExtendedAssociation.h:25
VertexFwd.h
TransientVertex.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
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
HLTDisplacedmumuVtxProducer::matchToPrevious_
const bool matchToPrevious_
Definition: HLTDisplacedmumuVtxProducer.h:44
TriggerRefsCollections.h
HLTDisplacedmumuVtxProducer.h
trigger
Definition: HLTPrescaleTableCond.h:8
HLTDisplacedmumuVtxProducer::minInvMass_
const double minInvMass_
Definition: HLTDisplacedmumuVtxProducer.h:48
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:183
HLTDisplacedmumuVtxProducer::minPt_
const double minPt_
Definition: HLTDisplacedmumuVtxProducer.h:46
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22
HLTDisplacedmumuVtxProducer::previousCandToken_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
Definition: HLTDisplacedmumuVtxProducer.h:43