CMS 3D CMS Logo

HLTmumutkVtxProducer.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <cmath>
3 #include <vector>
4 
8 
21 
22 #include "HLTmumutkVtxProducer.h"
24 
25 using namespace edm;
26 using namespace reco;
27 using namespace std;
28 using namespace trigger;
29 
30 // ----------------------------------------------------------------------
32  : muCandTag_(iConfig.getParameter<edm::InputTag>("MuCand")),
33  muCandToken_(consumes<reco::RecoChargedCandidateCollection>(muCandTag_)),
34  trkCandTag_(iConfig.getParameter<edm::InputTag>("TrackCand")),
35  trkCandToken_(consumes<reco::RecoChargedCandidateCollection>(trkCandTag_)),
36  previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
37  previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
38  mfName_(iConfig.getParameter<std::string>("SimpleMagneticField")),
39  thirdTrackMass_(iConfig.getParameter<double>("ThirdTrackMass")),
40  maxEta_(iConfig.getParameter<double>("MaxEta")),
41  minPt_(iConfig.getParameter<double>("MinPt")),
42  minInvMass_(iConfig.getParameter<double>("MinInvMass")),
43  maxInvMass_(iConfig.getParameter<double>("MaxInvMass")),
44  minD0Significance_(iConfig.getParameter<double>("MinD0Significance")),
45  overlapDR_(iConfig.getParameter<double>("OverlapDR")),
46  beamSpotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
47  beamSpotToken_(consumes<reco::BeamSpot>(beamSpotTag_)) {
48  produces<VertexCollection>();
49 }
50 
51 // ----------------------------------------------------------------------
53 
56  desc.add<edm::InputTag>("MuCand", edm::InputTag("hltMuTracks"));
57  desc.add<edm::InputTag>("TrackCand", edm::InputTag("hltMumukAllConeTracks"));
58  desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag("hltDisplacedmumuFilterDoubleMu4Jpsi"));
59  desc.add<std::string>("SimpleMagneticField", "");
60  desc.add<double>("ThirdTrackMass", 0.493677);
61  desc.add<double>("MaxEta", 2.5);
62  desc.add<double>("MinPt", 3.0);
63  desc.add<double>("MinInvMass", 0.0);
64  desc.add<double>("MaxInvMass", 99999.);
65  desc.add<double>("MinD0Significance", 0.0);
66  desc.add<double>("OverlapDR", 1.44e-4);
67  desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltOfflineBeamSpot"));
68  descriptions.add("HLTmumutkVtxProducer", desc);
69 }
70 
71 // ----------------------------------------------------------------------
73  const double MuMass(0.106);
74  const double MuMass2(MuMass * MuMass);
75  const double thirdTrackMass2(thirdTrackMass_ * thirdTrackMass_);
76 
77  // get hold of muon trks
79  iEvent.getByToken(muCandToken_, mucands);
80 
81  //get the transient track builder:
83  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
84 
85  //get the beamspot position
86  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
87  iEvent.getByToken(beamSpotToken_, recoBeamSpotHandle);
88 
89  //get the b field
90  ESHandle<MagneticField> bFieldHandle;
91  iSetup.get<IdealMagneticFieldRecord>().get(mfName_, bFieldHandle);
92  const MagneticField* magField = bFieldHandle.product();
93  TSCBLBuilderNoMaterial blsBuilder;
94 
95  // get track candidates around displaced muons
97  iEvent.getByToken(trkCandToken_, trkcands);
98 
99  unique_ptr<VertexCollection> vertexCollection(new VertexCollection());
100 
101  // Ref to Candidate object to be recorded in filter object
105 
106  double e1, e2, e3;
108 
109  if (mucands->size() < 2)
110  return;
111  if (trkcands->empty())
112  return;
113 
114  RecoChargedCandidateCollection::const_iterator mucand1;
115  RecoChargedCandidateCollection::const_iterator mucand2;
116  RecoChargedCandidateCollection::const_iterator trkcand;
117 
118  // get the objects passing the previous filter
120  iEvent.getByToken(previousCandToken_, previousCands);
121 
122  vector<RecoChargedCandidateRef> vPrevCands;
123  previousCands->getObjects(TriggerMuon, vPrevCands);
124 
125  for (mucand1 = mucands->begin(); mucand1 != mucands->end(); ++mucand1) {
126  TrackRef trk1 = mucand1->get<TrackRef>();
127  LogDebug("HLTmumutkVtxProducer") << " 1st muon: q*pt= " << trk1->charge() * trk1->pt() << ", eta= " << trk1->eta()
128  << ", hits= " << trk1->numberOfValidHits();
129 
130  //first check if this muon passed the previous filter
131  if (!checkPreviousCand(trk1, vPrevCands))
132  continue;
133 
134  // eta and pt cut
135  if (fabs(trk1->eta()) > maxEta_)
136  continue;
137  if (trk1->pt() < minPt_)
138  continue;
139 
140  mucand2 = mucand1;
141  ++mucand2;
142  for (; mucand2 != mucands->end(); mucand2++) {
143  TrackRef trk2 = mucand2->get<TrackRef>();
144 
145  LogDebug("HLTDisplacedMumukFilter") << " 2nd muon: q*pt= " << trk2->charge() * trk2->pt()
146  << ", eta= " << trk2->eta() << ", hits= " << trk2->numberOfValidHits();
147 
148  //first check if this muon passed the previous filter
149  if (!checkPreviousCand(trk2, vPrevCands))
150  continue;
151  // eta and pt cut
152  if (fabs(trk2->eta()) > maxEta_)
153  continue;
154  if (trk2->pt() < minPt_)
155  continue;
156 
157  //loop on track collection
158  for (trkcand = trkcands->begin(); trkcand != trkcands->end(); ++trkcand) {
159  TrackRef trk3 = trkcand->get<TrackRef>();
160  if (overlap(trk1, trk3))
161  continue;
162  if (overlap(trk2, trk3))
163  continue;
164 
165  LogDebug("HLTDisplacedMumukFilter") << " 3rd track: q*pt= " << trk3->charge() * trk3->pt()
166  << ", eta= " << trk3->eta() << ", hits= " << trk3->numberOfValidHits();
167 
168  // eta and pt cut
169  if (fabs(trk3->eta()) > maxEta_)
170  continue;
171  if (trk3->pt() < minPt_)
172  continue;
173 
174  // Combined system
175  e1 = sqrt(trk1->momentum().Mag2() + MuMass2);
176  e2 = sqrt(trk2->momentum().Mag2() + MuMass2);
177  e3 = sqrt(trk3->momentum().Mag2() + thirdTrackMass2);
178 
179  p1 = Particle::LorentzVector(trk1->px(), trk1->py(), trk1->pz(), e1);
180  p2 = Particle::LorentzVector(trk2->px(), trk2->py(), trk2->pz(), e2);
181  p3 = Particle::LorentzVector(trk3->px(), trk3->py(), trk3->pz(), e3);
182 
183  p = p1 + p2 + p3;
184 
185  //invariant mass cut
186  double invmass = abs(p.mass());
187  LogDebug("HLTDisplacedMumukFilter") << " Invmass= " << invmass;
188  if (invmass < minInvMass_)
189  continue;
190  if (invmass > maxInvMass_)
191  continue;
192 
193  // do the vertex fit
194  vector<TransientTrack> t_tks;
195  t_tks.push_back((*theB).build(&trk1));
196  t_tks.push_back((*theB).build(&trk2));
197  t_tks.push_back((*theB).build(&trk3));
198  if (t_tks.size() != 3)
199  continue;
200 
201  FreeTrajectoryState InitialFTS = initialFreeState(*trk3, magField);
202  TrajectoryStateClosestToBeamLine tscb(blsBuilder(InitialFTS, *recoBeamSpotHandle));
203  double d0sig = tscb.transverseImpactParameter().significance();
204  if (d0sig < minD0Significance_)
205  continue;
206 
207  KalmanVertexFitter kvf;
208  TransientVertex tv = kvf.vertex(t_tks);
209  if (!tv.isValid())
210  continue;
211  Vertex vertex = tv;
212 
213  // put vertex in the event
214  vertexCollection->push_back(vertex);
215  }
216  }
217  }
219 }
220 
223  GlobalPoint gpos(pos);
224  Basic3DVector<float> mom(tk.momentum());
225  GlobalVector gmom(mom);
226  GlobalTrajectoryParameters par(gpos, gmom, tk.charge(), field);
228  return FreeTrajectoryState(par, err);
229 }
230 
231 bool HLTmumutkVtxProducer::overlap(const TrackRef& trackref1, const TrackRef& trackref2) {
232  if (deltaR(trackref1->eta(), trackref1->phi(), trackref2->eta(), trackref2->phi()) < overlapDR_)
233  return true;
234  return false;
235 }
236 
238  const vector<RecoChargedCandidateRef>& refVect) const {
239  bool ok = false;
240  for (auto& i : refVect) {
241  if (i->get<TrackRef>() == trackref) {
242  ok = true;
243  break;
244  }
245  }
246  return ok;
247 }
Vector3DBase
Definition: Vector3DBase.h:8
ConfigurationDescriptions.h
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
TSCBLBuilderNoMaterial.h
TrajectoryStateClosestToBeamLine
Definition: TrajectoryStateClosestToBeamLine.h:15
trigger::TriggerFilterObjectWithRefs
Definition: TriggerFilterObjectWithRefs.h:35
KalmanVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition: KalmanVertexFitter.h:49
HLTmumutkVtxProducer.h
mps_fire.i
i
Definition: mps_fire.py:428
FreeTrajectoryState.h
GlobalTrajectoryParameters.h
MessageLogger.h
HLTmumutkVtxProducer::mfName_
const std::string mfName_
Definition: HLTmumutkVtxProducer.h:62
KalmanVertexFitter.h
ESHandle.h
HLTmumutkVtxProducer::initialFreeState
static FreeTrajectoryState initialFreeState(const reco::Track &, const MagneticField *)
Definition: HLTmumutkVtxProducer.cc:221
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
pos
Definition: PixelAliasList.h:18
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
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
HLTmumutkVtxProducer::minPt_
const double minPt_
Definition: HLTmumutkVtxProducer.h:66
edm::Ref
Definition: AssociativeIterator.h:58
HLTmumutkVtxProducer::minD0Significance_
const double minD0Significance_
Definition: HLTmumutkVtxProducer.h:69
IdealMagneticFieldRecord
Definition: IdealMagneticFieldRecord.h:11
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
CurvilinearTrajectoryError
Definition: CurvilinearTrajectoryError.h:27
BeamSpot.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
reco::BeamSpot
Definition: BeamSpot.h:21
TransientTrackRecord
Definition: TransientTrackRecord.h:11
trigger::TriggerMuon
Definition: TriggerTypeDefs.h:68
reco::Track
Definition: Track.h:27
edm::ESHandle< TransientTrackBuilder >
p2
double p2[4]
Definition: TauolaWrapper.h:90
reco::TrackBase::charge
int charge() const
track electric charge
Definition: TrackBase.h:596
GlobalTrajectoryParameters
Definition: GlobalTrajectoryParameters.h:15
TrajectoryStateClosestToBeamLine::transverseImpactParameter
Measurement1D transverseImpactParameter() const
Definition: TrajectoryStateClosestToBeamLine.cc:3
reco::TrackBase::covariance
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:716
Point3DBase< float, GlobalTag >
HLTmumutkVtxProducer::maxEta_
const double maxEta_
Definition: HLTmumutkVtxProducer.h:65
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
Measurement1D::significance
double significance() const
Definition: Measurement1D.h:29
bphysicsOniaDQM_cfi.vertex
vertex
Definition: bphysicsOniaDQM_cfi.py:7
TransientTrackBuilder.h
HLTmumutkVtxProducer::~HLTmumutkVtxProducer
~HLTmumutkVtxProducer() override
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
edm::ParameterSet
Definition: ParameterSet.h:47
deltaR.h
StorageManager_cfg.e1
e1
Definition: StorageManager_cfg.py:16
HLTmumutkVtxProducer::minInvMass_
const double minInvMass_
Definition: HLTmumutkVtxProducer.h:67
iEvent
int iEvent
Definition: GenABIO.cc:224
HLTmumutkVtxProducer::checkPreviousCand
bool checkPreviousCand(const reco::TrackRef &trackref, const std::vector< reco::RecoChargedCandidateRef > &ref2) const
Definition: HLTmumutkVtxProducer.cc:237
p1
double p1[4]
Definition: TauolaWrapper.h:89
HLTmumutkVtxProducer::previousCandToken_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
Definition: HLTmumutkVtxProducer.h:60
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
TransientTrackRecord.h
submitPVResolutionJobs.err
err
Definition: submitPVResolutionJobs.py:85
get
#define get
reco::TrackBase::vertex
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead
Definition: TrackBase.h:676
HLTmumutkVtxProducer::beamSpotToken_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
Definition: HLTmumutkVtxProducer.h:73
SiPixelPhase1Clusters_cfi.e3
e3
Definition: SiPixelPhase1Clusters_cfi.py:9
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
FreeTrajectoryState
Definition: FreeTrajectoryState.h:27
HLTmumutkVtxProducer::thirdTrackMass_
const double thirdTrackMass_
Definition: HLTmumutkVtxProducer.h:64
RecoChargedCandidate.h
TSCBLBuilderNoMaterial
Definition: TSCBLBuilderNoMaterial.h:13
spclusmultinvestigator_cfi.vertexCollection
vertexCollection
Definition: spclusmultinvestigator_cfi.py:4
HLTmumutkVtxProducer::maxInvMass_
const double maxInvMass_
Definition: HLTmumutkVtxProducer.h:68
HLTmumutkVtxProducer::HLTmumutkVtxProducer
HLTmumutkVtxProducer(const edm::ParameterSet &)
Definition: HLTmumutkVtxProducer.cc:31
HLTmumutkVtxProducer::trkCandToken_
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > trkCandToken_
Definition: HLTmumutkVtxProducer.h:58
p3
double p3[4]
Definition: TauolaWrapper.h:91
trigger
Definition: HLTPrescaleTableCond.h:8
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::TrackBase::momentum
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
HLTmumutkVtxProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTmumutkVtxProducer.cc:54
edm::Event
Definition: Event.h:73
MagneticField
Definition: MagneticField.h:19
HLTmumutkVtxProducer::overlapDR_
const double overlapDR_
Definition: HLTmumutkVtxProducer.h:70
HLTmumutkVtxProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: HLTmumutkVtxProducer.cc:72
HLTmumutkVtxProducer::overlap
bool overlap(const reco::TrackRef &trackref1, const reco::TrackRef &trackref2)
Definition: HLTmumutkVtxProducer.cc:231
edm::InputTag
Definition: InputTag.h:15
Basic3DVector< float >
HLTmumutkVtxProducer::muCandToken_
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > muCandToken_
Definition: HLTmumutkVtxProducer.h:56
reco::Vertex
Definition: Vertex.h:35
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22