CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
HLTmumutktkVtxProducer Class Reference

#include <HLTmumutktkVtxProducer.h>

Inheritance diagram for HLTmumutktkVtxProducer:
edm::stream::EDProducer<>

Public Member Functions

 HLTmumutktkVtxProducer (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~HLTmumutktkVtxProducer () override=default
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

bool checkPreviousCand (const reco::TrackRef &trackref, const std::vector< reco::RecoChargedCandidateRef > &ref2) const
 
bool overlap (const reco::TrackRef &trackref1, const reco::TrackRef &trackref2)
 

Static Private Member Functions

static FreeTrajectoryState initialFreeState (const reco::Track &, const MagneticField *)
 

Private Attributes

const edm::InputTag beamSpotTag_
 
const edm::EDGetTokenT< reco::BeamSpotbeamSpotToken_
 
const double fourthTrackMass_
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordidealMagneticFieldRecordToken_
 
const double maxEta_
 
const double maxInvMass_
 
const double maxTrkTrkMass_
 
const std::string mfName_
 
const double minD0Significance_
 
const double minInvMass_
 
const double minPt_
 
const double minTrkTrkMass_
 
const edm::InputTag muCandTag_
 
const edm::EDGetTokenT< reco::RecoChargedCandidateCollectionmuCandToken_
 
const bool oppositeSign_
 
const double overlapDR2_
 
const edm::InputTag previousCandTag_
 
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefspreviousCandToken_
 
const double thirdTrackMass_
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordtransientTrackRecordToken_
 
const edm::InputTag trkCandTag_
 
const edm::EDGetTokenT< reco::RecoChargedCandidateCollectiontrkCandToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 42 of file HLTmumutktkVtxProducer.h.

Constructor & Destructor Documentation

◆ HLTmumutktkVtxProducer()

HLTmumutktkVtxProducer::HLTmumutktkVtxProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 26 of file HLTmumutktkVtxProducer.cc.

27  : transientTrackRecordToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
28  muCandTag_(iConfig.getParameter<edm::InputTag>("MuCand")),
29  muCandToken_(consumes<reco::RecoChargedCandidateCollection>(muCandTag_)),
30  trkCandTag_(iConfig.getParameter<edm::InputTag>("TrackCand")),
31  trkCandToken_(consumes<reco::RecoChargedCandidateCollection>(trkCandTag_)),
32  previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
33  previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
34  mfName_(iConfig.getParameter<std::string>("SimpleMagneticField")),
36  thirdTrackMass_(iConfig.getParameter<double>("ThirdTrackMass")),
37  fourthTrackMass_(iConfig.getParameter<double>("FourthTrackMass")),
38  maxEta_(iConfig.getParameter<double>("MaxEta")),
39  minPt_(iConfig.getParameter<double>("MinPt")),
40  minInvMass_(iConfig.getParameter<double>("MinInvMass")),
41  maxInvMass_(iConfig.getParameter<double>("MaxInvMass")),
42  minTrkTrkMass_(iConfig.getParameter<double>("MinTrkTrkMass")),
43  maxTrkTrkMass_(iConfig.getParameter<double>("MaxTrkTrkMass")),
44  minD0Significance_(iConfig.getParameter<double>("MinD0Significance")),
45  oppositeSign_(iConfig.getParameter<bool>("OppositeSign")),
46  // minimum delta-R^2 threshold (with sign) for non-overlapping tracks
47  overlapDR2_(iConfig.getParameter<double>("OverlapDR") * std::abs(iConfig.getParameter<double>("OverlapDR"))),
48  beamSpotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
49  beamSpotToken_(consumes<reco::BeamSpot>(beamSpotTag_)) {
50  produces<VertexCollection>();
51 }
const edm::InputTag previousCandTag_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::InputTag muCandTag_
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > trkCandToken_
const edm::InputTag beamSpotTag_
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > muCandToken_
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackRecordToken_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const edm::InputTag trkCandTag_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > idealMagneticFieldRecordToken_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_

◆ ~HLTmumutktkVtxProducer()

HLTmumutktkVtxProducer::~HLTmumutktkVtxProducer ( )
overridedefault

Member Function Documentation

◆ checkPreviousCand()

bool HLTmumutktkVtxProducer::checkPreviousCand ( const reco::TrackRef trackref,
const std::vector< reco::RecoChargedCandidateRef > &  ref2 
) const
private

Definition at line 280 of file HLTmumutktkVtxProducer.cc.

References mps_fire::i, and convertSQLiteXML::ok.

Referenced by produce().

281  {
282  bool ok = false;
283  for (auto& i : refVect) {
284  if (i->get<TrackRef>() == trackref) {
285  ok = true;
286  break;
287  }
288  }
289  return ok;
290 }

◆ fillDescriptions()

void HLTmumutktkVtxProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 53 of file HLTmumutktkVtxProducer.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

53  {
55  desc.add<edm::InputTag>("MuCand", edm::InputTag("hltMuTracks"));
56  desc.add<edm::InputTag>("TrackCand", edm::InputTag("hltMumukAllConeTracks"));
57  desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag("hltDisplacedmumuFilterDoubleMu4Jpsi"));
58  desc.add<std::string>("SimpleMagneticField", "");
59  desc.add<double>("ThirdTrackMass", 0.493677);
60  desc.add<double>("FourthTrackMass", 0.493677);
61  desc.add<double>("MaxEta", 2.5);
62  desc.add<double>("MinPt", 0.0);
63  desc.add<double>("MinInvMass", 0.0);
64  desc.add<double>("MaxInvMass", 99999.);
65  desc.add<double>("MinTrkTrkMass", 0.0);
66  desc.add<double>("MaxTrkTrkMass", 99999.);
67  desc.add<double>("MinD0Significance", 0.0);
68  desc.add<bool>("OppositeSign", false);
69  desc.add<double>("OverlapDR", 0.001);
70  desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltOfflineBeamSpot"));
71  descriptions.add("HLTmumutktkVtxProducer", desc);
72 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ initialFreeState()

FreeTrajectoryState HLTmumutktkVtxProducer::initialFreeState ( const reco::Track tk,
const MagneticField field 
)
staticprivate

Definition at line 266 of file HLTmumutktkVtxProducer.cc.

References reco::TrackBase::charge(), reco::TrackBase::covariance(), submitPVResolutionJobs::err, reco::TrackBase::momentum(), and reco::TrackBase::vertex().

Referenced by produce().

266  {
268  GlobalPoint gpos(pos);
269  Basic3DVector<float> mom(tk.momentum());
270  GlobalVector gmom(mom);
271  GlobalTrajectoryParameters par(gpos, gmom, tk.charge(), field);
273  return FreeTrajectoryState(par, err);
274 }
int charge() const
track electric charge
Definition: TrackBase.h:596
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:716
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:676
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664

◆ overlap()

bool HLTmumutktkVtxProducer::overlap ( const reco::TrackRef trackref1,
const reco::TrackRef trackref2 
)
private

Definition at line 276 of file HLTmumutktkVtxProducer.cc.

References reco::deltaR2(), and overlapDR2_.

Referenced by produce().

276  {
277  return (reco::deltaR2(trackref1->eta(), trackref1->phi(), trackref2->eta(), trackref2->phi()) < overlapDR2_);
278 }
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16

◆ produce()

void HLTmumutktkVtxProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 74 of file HLTmumutktkVtxProducer.cc.

References beamSpotToken_, checkPreviousCand(), StorageManager_cfg::e1, fourthTrackMass_, edm::Ref< C, T, F >::get(), edm::EventSetup::getHandle(), trigger::TriggerRefsCollections::getObjects(), idealMagneticFieldRecordToken_, iEvent, initialFreeState(), TransientVertex::isValid(), LogDebug, maxEta_, maxInvMass_, maxTrkTrkMass_, minD0Significance_, minInvMass_, minPt_, minTrkTrkMass_, eostools::move(), muCandToken_, oppositeSign_, overlap(), AlCaHLTBitMon_ParallelJobs::p, LaserDQM_cfg::p1, SiStripOfflineCRack_cfg::p2, previousCandToken_, Measurement1D::significance(), mathSSE::sqrt(), thirdTrackMass_, transientTrackRecordToken_, TrajectoryStateClosestToBeamLine::transverseImpactParameter(), trigger::TriggerMuon, trkCandToken_, bphysicsOniaDQM_cfi::vertex, KalmanVertexFitter::vertex(), and spclusmultinvestigator_cfi::vertexCollection.

74  {
75  const double MuMass(0.106);
76  const double MuMass2(MuMass * MuMass);
77  const double thirdTrackMass2(thirdTrackMass_ * thirdTrackMass_);
78  const double fourthTrackMass2(fourthTrackMass_ * fourthTrackMass_);
79 
80  // get hold of muon trks
82  iEvent.getByToken(muCandToken_, mucands);
83 
84  // get the transient track builder
85  auto const& theB = iSetup.getHandle(transientTrackRecordToken_);
86 
87  // get the beamspot position
88  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
89  iEvent.getByToken(beamSpotToken_, recoBeamSpotHandle);
90 
91  // get the b field
92  auto const& bFieldHandle = iSetup.getHandle(idealMagneticFieldRecordToken_);
93  const MagneticField* magField = bFieldHandle.product();
94  TSCBLBuilderNoMaterial blsBuilder;
95 
96  // get track candidates around displaced muons
98  iEvent.getByToken(trkCandToken_, trkcands);
99 
100  unique_ptr<VertexCollection> vertexCollection(new VertexCollection());
101 
102  // Ref to Candidate object to be recorded in filter object
105  RecoChargedCandidateRef refTrk1;
106  RecoChargedCandidateRef refTrk2;
107 
108  double e1, e2, e3_m3, e3_m4, e4_m3, e4_m4;
109  Particle::LorentzVector p, pBar, p1, p2, p3_m3, p3_m4, p4_m3, p4_m4, p_m3m4, p_m4m3;
110 
111  if (mucands->size() < 2)
112  return;
113  if (trkcands->size() < 2)
114  return;
115 
116  RecoChargedCandidateCollection::const_iterator mucand1;
117  RecoChargedCandidateCollection::const_iterator mucand2;
118  RecoChargedCandidateCollection::const_iterator trkcand1;
119  RecoChargedCandidateCollection::const_iterator trkcand2;
120 
121  // get the objects passing the previous filter
123  iEvent.getByToken(previousCandToken_, previousCands);
124 
125  vector<RecoChargedCandidateRef> vPrevCands;
126  previousCands->getObjects(TriggerMuon, vPrevCands);
127 
128  for (mucand1 = mucands->begin(); mucand1 != mucands->end(); ++mucand1) {
129  TrackRef trk1 = mucand1->get<TrackRef>();
130  LogDebug("HLTmumutktkVtxProducer") << " 1st muon: q*pt= " << trk1->charge() * trk1->pt() << ", eta= " << trk1->eta()
131  << ", hits= " << trk1->numberOfValidHits();
132 
133  //first check if this muon passed the previous filter
134  if (!checkPreviousCand(trk1, vPrevCands))
135  continue;
136  // eta and pt cut
137  if (fabs(trk1->eta()) > maxEta_)
138  continue;
139  if (trk1->pt() < minPt_)
140  continue;
141 
142  mucand2 = mucand1;
143  ++mucand2;
144  for (; mucand2 != mucands->end(); mucand2++) {
145  TrackRef trk2 = mucand2->get<TrackRef>();
146  if (overlap(trk1, trk2))
147  continue;
148 
149  LogDebug("HLTDisplacedMumukFilter") << " 2nd muon: q*pt= " << trk2->charge() * trk2->pt()
150  << ", eta= " << trk2->eta() << ", hits= " << trk2->numberOfValidHits();
151 
152  //first check if this muon passed the previous filter
153  if (!checkPreviousCand(trk2, vPrevCands))
154  continue;
155  // eta and pt cut
156  if (fabs(trk2->eta()) > maxEta_)
157  continue;
158  if (trk2->pt() < minPt_)
159  continue;
160 
161  //loop on track collection - trk1
162  for (trkcand1 = trkcands->begin(); trkcand1 != trkcands->end(); ++trkcand1) {
163  TrackRef trk3 = trkcand1->get<TrackRef>();
164 
165  if (overlap(trk1, trk3))
166  continue;
167  if (overlap(trk2, trk3))
168  continue;
169 
170  LogDebug("HLTDisplacedMumukFilter") << " 3rd track: q*pt= " << trk3->charge() * trk3->pt()
171  << ", eta= " << trk3->eta() << ", hits= " << trk3->numberOfValidHits();
172 
173  // eta and pt cut
174  if (fabs(trk3->eta()) > maxEta_)
175  continue;
176  if (trk3->pt() < minPt_)
177  continue;
178 
179  FreeTrajectoryState InitialFTS_Trk3 = initialFreeState(*trk3, magField);
180  TrajectoryStateClosestToBeamLine tscb_Trk3(blsBuilder(InitialFTS_Trk3, *recoBeamSpotHandle));
181  double d0sigTrk3 = tscb_Trk3.transverseImpactParameter().significance();
182  if (d0sigTrk3 < minD0Significance_)
183  continue;
184 
185  //loop on track collection - trk2
186  for (trkcand2 = trkcands->begin(); trkcand2 != trkcands->end(); ++trkcand2) {
187  TrackRef trk4 = trkcand2->get<TrackRef>();
188 
189  if (oppositeSign_) {
190  if (trk3->charge() * trk4->charge() != -1)
191  continue;
192  }
193  if (overlap(trk1, trk4))
194  continue;
195  if (overlap(trk2, trk4))
196  continue;
197  if (overlap(trk3, trk4))
198  continue;
199 
200  LogDebug("HLTDisplacedMumukFilter") << " 4th track: q*pt= " << trk4->charge() * trk4->pt()
201  << ", eta= " << trk4->eta() << ", hits= " << trk4->numberOfValidHits();
202 
203  // eta and pt cut
204  if (fabs(trk4->eta()) > maxEta_)
205  continue;
206  if (trk4->pt() < minPt_)
207  continue;
208 
209  FreeTrajectoryState InitialFTS_Trk4 = initialFreeState(*trk4, magField);
210  TrajectoryStateClosestToBeamLine tscb_Trk4(blsBuilder(InitialFTS_Trk4, *recoBeamSpotHandle));
211  double d0sigTrk4 = tscb_Trk4.transverseImpactParameter().significance();
212  if (d0sigTrk4 < minD0Significance_)
213  continue;
214 
215  // Combined system
216  e1 = sqrt(trk1->momentum().Mag2() + MuMass2);
217  e2 = sqrt(trk2->momentum().Mag2() + MuMass2);
218  e3_m3 = sqrt(trk3->momentum().Mag2() + thirdTrackMass2);
219  e3_m4 = sqrt(trk3->momentum().Mag2() + fourthTrackMass2);
220  e4_m3 = sqrt(trk4->momentum().Mag2() + thirdTrackMass2);
221  e4_m4 = sqrt(trk4->momentum().Mag2() + fourthTrackMass2);
222 
223  p1 = Particle::LorentzVector(trk1->px(), trk1->py(), trk1->pz(), e1);
224  p2 = Particle::LorentzVector(trk2->px(), trk2->py(), trk2->pz(), e2);
225  p3_m3 = Particle::LorentzVector(trk3->px(), trk3->py(), trk3->pz(), e3_m3);
226  p3_m4 = Particle::LorentzVector(trk3->px(), trk3->py(), trk3->pz(), e3_m4);
227  p4_m3 = Particle::LorentzVector(trk4->px(), trk4->py(), trk4->pz(), e4_m3);
228  p4_m4 = Particle::LorentzVector(trk4->px(), trk4->py(), trk4->pz(), e4_m4);
229 
230  p = p1 + p2 + p3_m3 + p4_m4;
231  pBar = p1 + p2 + p3_m4 + p4_m3;
232  p_m3m4 = p3_m3 + p4_m4;
233  p_m4m3 = p3_m4 + p4_m3;
234 
235  //invariant mass cut
236  if (!((p_m3m4.mass() > minTrkTrkMass_ && p_m3m4.mass() < maxTrkTrkMass_) ||
237  (p_m4m3.mass() > minTrkTrkMass_ && p_m4m3.mass() < maxTrkTrkMass_)))
238  continue;
239  if (!((p.mass() > minInvMass_ && p.mass() < maxInvMass_) ||
240  (pBar.mass() > minInvMass_ && pBar.mass() < maxInvMass_)))
241  continue;
242 
243  // do the vertex fit
244  vector<TransientTrack> t_tks;
245  t_tks.push_back((*theB).build(&trk1));
246  t_tks.push_back((*theB).build(&trk2));
247  t_tks.push_back((*theB).build(&trk3));
248  t_tks.push_back((*theB).build(&trk4));
249  if (t_tks.size() != 4)
250  continue;
251 
252  KalmanVertexFitter kvf;
253  TransientVertex tv = kvf.vertex(t_tks);
254  if (!tv.isValid())
255  continue;
256  Vertex vertex = tv;
257 
258  vertexCollection->push_back(vertex);
259  }
260  }
261  }
262  }
264 }
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool checkPreviousCand(const reco::TrackRef &trackref, const std::vector< reco::RecoChargedCandidateRef > &ref2) const
bool overlap(const reco::TrackRef &trackref1, const reco::TrackRef &trackref2)
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
static FreeTrajectoryState initialFreeState(const reco::Track &, const MagneticField *)
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > trkCandToken_
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
const edm::EDGetTokenT< reco::RecoChargedCandidateCollection > muCandToken_
int iEvent
Definition: GenABIO.cc:224
bool isValid() const
T sqrt(T t)
Definition: SSEVec.h:23
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transientTrackRecordToken_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > idealMagneticFieldRecordToken_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > previousCandToken_
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
def move(src, dest)
Definition: eostools.py:511
math::PtEtaPhiELorentzVectorF LorentzVector
#define LogDebug(id)

Member Data Documentation

◆ beamSpotTag_

const edm::InputTag HLTmumutktkVtxProducer::beamSpotTag_
private

Definition at line 79 of file HLTmumutktkVtxProducer.h.

◆ beamSpotToken_

const edm::EDGetTokenT<reco::BeamSpot> HLTmumutktkVtxProducer::beamSpotToken_
private

Definition at line 80 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ fourthTrackMass_

const double HLTmumutktkVtxProducer::fourthTrackMass_
private

Definition at line 69 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ idealMagneticFieldRecordToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> HLTmumutktkVtxProducer::idealMagneticFieldRecordToken_
private

Definition at line 66 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ maxEta_

const double HLTmumutktkVtxProducer::maxEta_
private

Definition at line 70 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ maxInvMass_

const double HLTmumutktkVtxProducer::maxInvMass_
private

Definition at line 73 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ maxTrkTrkMass_

const double HLTmumutktkVtxProducer::maxTrkTrkMass_
private

Definition at line 75 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ mfName_

const std::string HLTmumutktkVtxProducer::mfName_
private

Definition at line 65 of file HLTmumutktkVtxProducer.h.

◆ minD0Significance_

const double HLTmumutktkVtxProducer::minD0Significance_
private

Definition at line 76 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ minInvMass_

const double HLTmumutktkVtxProducer::minInvMass_
private

Definition at line 72 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ minPt_

const double HLTmumutktkVtxProducer::minPt_
private

Definition at line 71 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ minTrkTrkMass_

const double HLTmumutktkVtxProducer::minTrkTrkMass_
private

Definition at line 74 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ muCandTag_

const edm::InputTag HLTmumutktkVtxProducer::muCandTag_
private

Definition at line 58 of file HLTmumutktkVtxProducer.h.

◆ muCandToken_

const edm::EDGetTokenT<reco::RecoChargedCandidateCollection> HLTmumutktkVtxProducer::muCandToken_
private

Definition at line 59 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ oppositeSign_

const bool HLTmumutktkVtxProducer::oppositeSign_
private

Definition at line 77 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ overlapDR2_

const double HLTmumutktkVtxProducer::overlapDR2_
private

Definition at line 78 of file HLTmumutktkVtxProducer.h.

Referenced by overlap().

◆ previousCandTag_

const edm::InputTag HLTmumutktkVtxProducer::previousCandTag_
private

Definition at line 62 of file HLTmumutktkVtxProducer.h.

◆ previousCandToken_

const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> HLTmumutktkVtxProducer::previousCandToken_
private

Definition at line 63 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ thirdTrackMass_

const double HLTmumutktkVtxProducer::thirdTrackMass_
private

Definition at line 68 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ transientTrackRecordToken_

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> HLTmumutktkVtxProducer::transientTrackRecordToken_
private

Definition at line 56 of file HLTmumutktkVtxProducer.h.

Referenced by produce().

◆ trkCandTag_

const edm::InputTag HLTmumutktkVtxProducer::trkCandTag_
private

Definition at line 60 of file HLTmumutktkVtxProducer.h.

◆ trkCandToken_

const edm::EDGetTokenT<reco::RecoChargedCandidateCollection> HLTmumutktkVtxProducer::trkCandToken_
private

Definition at line 61 of file HLTmumutktkVtxProducer.h.

Referenced by produce().