CMS 3D CMS Logo

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

#include <TransientTrackBuilder.h>

Public Member Functions

reco::TransientTrack build (const reco::Track *p) const
 
reco::TransientTrack build (const reco::Track &p) const
 
reco::TransientTrack build (const reco::GsfTrack *p) const
 
reco::TransientTrack build (const reco::GsfTrack &p) const
 
reco::TransientTrack build (const reco::TrackRef *p) const
 
reco::TransientTrack build (const reco::TrackRef &p) const
 
reco::TransientTrack build (const reco::GsfTrackRef *p) const
 
reco::TransientTrack build (const reco::GsfTrackRef &p) const
 
reco::TransientTrack build (const reco::CandidatePtr *p) const
 
reco::TransientTrack build (const reco::CandidatePtr &p) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::TrackCollection > &trkColl) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::GsfTrackCollection > &trkColl) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< edm::View< reco::Track > > &trkColl) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::TrackCollection > &trkColl, const edm::ValueMap< float > &trackTimes, const edm::ValueMap< float > &trackTimeResos) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::GsfTrackCollection > &trkColl, const edm::ValueMap< float > &trackTimes, const edm::ValueMap< float > &trackTimeResos) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< edm::View< reco::Track > > &trkColl, const edm::ValueMap< float > &trackTimes, const edm::ValueMap< float > &trackTimeResos) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::TrackCollection > &trkColl, const reco::BeamSpot &beamSpot) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::GsfTrackCollection > &trkColl, const reco::BeamSpot &beamSpot) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< edm::View< reco::Track > > &trkColl, const reco::BeamSpot &beamSpot) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::TrackCollection > &trkColl, const reco::BeamSpot &beamSpot, const edm::ValueMap< float > &trackTimes, const edm::ValueMap< float > &trackTimeResos) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< reco::GsfTrackCollection > &trkColl, const reco::BeamSpot &beamSpot, const edm::ValueMap< float > &trackTimes, const edm::ValueMap< float > &trackTimeResos) const
 
std::vector< reco::TransientTrackbuild (const edm::Handle< edm::View< reco::Track > > &trkColl, const reco::BeamSpot &beamSpot, const edm::ValueMap< float > &trackTimes, const edm::ValueMap< float > &trackTimeResos) const
 
reco::TransientTrack build (const FreeTrajectoryState &fts) const
 
const MagneticFieldfield () const
 
const edm::ESHandle< GlobalTrackingGeometrytrackingGeometry () const
 
 TransientTrackBuilder (const MagneticField *field, const edm::ESHandle< GlobalTrackingGeometry > &trackingGeometry)
 

Static Public Attributes

static constexpr float defaultInvalidTrackTimeReso = 0.350f
 

Private Attributes

const MagneticFieldtheField
 
edm::ESHandle< GlobalTrackingGeometrytheTrackingGeometry
 

Detailed Description

Helper class to build TransientTrack from the persistent Track. This is obtained from the eventSetup, as given in the example in the test directory.

Definition at line 16 of file TransientTrackBuilder.h.

Constructor & Destructor Documentation

◆ TransientTrackBuilder()

TransientTrackBuilder::TransientTrackBuilder ( const MagneticField field,
const edm::ESHandle< GlobalTrackingGeometry > &  trackingGeometry 
)
inline

Definition at line 18 of file TransientTrackBuilder.h.

const MagneticField * theField
const MagneticField * field() const
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
const edm::ESHandle< GlobalTrackingGeometry > trackingGeometry() const

Member Function Documentation

◆ build() [1/23]

TransientTrack TransientTrackBuilder::build ( const reco::Track p) const

Definition at line 16 of file TransientTrackBuilder.cc.

References submitPVValidationJobs::t.

Referenced by MuonEnergyDepositAnalyzer::analyze(), DiMuonVertexMonitor::analyze(), tadqm::TrackAnalyzer::analyze(), SiPixelTrackResidualSource::analyze(), StandaloneTrackMonitor::analyze(), GEMEfficiencyAnalyzer::analyze(), DiMuonVertexValidation::analyze(), DiElectronVertexValidation::analyze(), PrimaryVertexValidation::analyze(), SplitVertexResolution::analyze(), HIPTwoBodyDecayAnalyzer::analyzeTrackCollection(), ConversionProducer::buildCollection(), ConvertedPhotonProducer::buildCollections(), btagbtvdeep::TrackInfoBuilder::buildTrackInfo(), tthelpers::buildTT(), pat::LeptonVertexSignificance::calculate(), PrimaryVertexAssignment::chargedHadronVertex(), reco::JetSignalVertexCompatibilityAlgo::convert(), DiMuonMassBiasMonitor::fillDecayHistograms(), tadqm::TrackAnalyzer::fillHistosForState(), DeepBoostedJetTagInfoProducer::fillParticleFeatures(), MuGEMMuonExtTableProducer::fillTable(), DisappearingMuonsSkimming::filter(), AlignmentTrackFromVertexSelector::findClosestVertex(), HIPTwoBodyDecayAnalyzer::fitDimuonVertex(), BPHDecayVertex::fTTracks(), EGammaMvaEleEstimator::IDIsoCombinedMvaValue(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), OniaVtxReProducer::makeVertices(), EGammaMvaEleEstimatorCSA14::mvaValue(), EGammaMvaEleEstimator::mvaValue(), reco::tau::RecoTauImpactParameterSignificancePlugin::operator()(), CheckHitPattern::operator()(), PFTrackProducer::produce(), MuonBeamspotConstraintValueMapProducer::produce(), Onia2MuMuPAT::produce(), ConvertedPhotonProducer::produce(), SoftPFMuonTagInfoProducer::produce(), BoostedDoubleSVProducer::produce(), ParticleTransformerAK4TagInfoProducer::produce(), pat::PATElectronProducer::produce(), DeepFlavourTagInfoProducer::produce(), pat::PATMuonProducer::produce(), TemplatedSecondaryVertexProducer< IPTI, VTX >::produce(), IPProducer< Container, Base, Helper >::produce(), CSCOverlapsAlignmentAlgorithm::run(), ConvBremPFTrackFinder::runConvBremFinder(), PrimaryVertexResolution::sortTracksByPt(), SoftLepton::tag(), TrackEfficiencyMonitor::testSTATracks(), and TrackEfficiencyMonitor::testTrackerTracks().

16  {
18 }
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

◆ build() [2/23]

TransientTrack TransientTrackBuilder::build ( const reco::Track p) const

Definition at line 20 of file TransientTrackBuilder.cc.

References submitPVValidationJobs::t.

20  {
22 }
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

◆ build() [3/23]

TransientTrack TransientTrackBuilder::build ( const reco::GsfTrack p) const

◆ build() [4/23]

TransientTrack TransientTrackBuilder::build ( const reco::GsfTrack p) const

◆ build() [5/23]

TransientTrack TransientTrackBuilder::build ( const reco::TrackRef p) const

Definition at line 45 of file TransientTrackBuilder.cc.

References submitPVValidationJobs::t.

45  {
47 }
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

◆ build() [6/23]

TransientTrack TransientTrackBuilder::build ( const reco::TrackRef p) const

Definition at line 49 of file TransientTrackBuilder.cc.

References submitPVValidationJobs::t.

49  {
51 }
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

◆ build() [7/23]

TransientTrack TransientTrackBuilder::build ( const reco::GsfTrackRef p) const

◆ build() [8/23]

TransientTrack TransientTrackBuilder::build ( const reco::GsfTrackRef p) const

◆ build() [9/23]

TransientTrack TransientTrackBuilder::build ( const reco::CandidatePtr p) const

Definition at line 32 of file TransientTrackBuilder.cc.

References edm::Ptr< T >::get(), reco::PFCandidate::isTimeValid(), submitPVValidationJobs::t, reco::PFCandidate::time(), pat::PackedCandidate::time(), reco::PFCandidate::timeError(), and pat::PackedCandidate::timeError().

32  {
33  reco::PFCandidatePtr tryPF(*t);
35  if (tryPF.get() != nullptr && tryPF->isTimeValid()) {
36  return TransientTrack(*t, tryPF->time(), tryPF->timeError(), theField, theTrackingGeometry);
37  } else if (tryPacked.get() != nullptr && tryPacked->timeError() > 0.f) {
38  return TransientTrack(*t, (double)tryPacked->time(), (double)tryPacked->timeError(), theField, theTrackingGeometry);
39  }
41 }
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

◆ build() [10/23]

TransientTrack TransientTrackBuilder::build ( const reco::CandidatePtr p) const

Definition at line 43 of file TransientTrackBuilder.cc.

References newFWLiteAna::build.

43 { return this->build(&t); }
reco::TransientTrack build(const reco::Track *p) const

◆ build() [11/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::TrackCollection > &  trkColl) const

Definition at line 61 of file TransientTrackBuilder.cc.

References mps_fire::i.

61  {
62  vector<TransientTrack> ttVect;
63  ttVect.reserve((*trkColl).size());
64  for (unsigned int i = 0; i < (*trkColl).size(); i++) {
65  ttVect.push_back(TransientTrack(TrackRef(trkColl, i), theField, theTrackingGeometry));
66  }
67  return ttVect;
68 }
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20

◆ build() [12/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::GsfTrackCollection > &  trkColl) const

Definition at line 70 of file TransientTrackBuilder.cc.

References mps_fire::i.

70  {
71  vector<TransientTrack> ttVect;
72  ttVect.reserve((*trkColl).size());
73  for (unsigned int i = 0; i < (*trkColl).size(); i++) {
74  ttVect.push_back(TransientTrack(new GsfTransientTrack(GsfTrackRef(trkColl, i), theField, theTrackingGeometry)));
75  }
76  return ttVect;
77 }
edm::Ref< GsfTrackCollection > GsfTrackRef
persistent reference to a GsfTrack
Definition: GsfTrackFwd.h:13
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

◆ build() [13/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< edm::View< reco::Track > > &  trkColl) const

Definition at line 79 of file TransientTrackBuilder.cc.

References mps_fire::i.

79  {
80  vector<TransientTrack> ttVect;
81  ttVect.reserve((*trkColl).size());
82  for (unsigned int i = 0; i < (*trkColl).size(); i++) {
83  const Track* trk = &(*trkColl)[i];
84  const GsfTrack* gsfTrack = dynamic_cast<const GsfTrack*>(trk);
85  if (gsfTrack) {
86  ttVect.push_back(TransientTrack(
87  new GsfTransientTrack(RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>(), theField, theTrackingGeometry)));
88  } else { // gsf
89  ttVect.push_back(TransientTrack(RefToBase<Track>(trkColl, i).castTo<TrackRef>(), theField, theTrackingGeometry));
90  }
91  }
92  return ttVect;
93 }
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

◆ build() [14/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::TrackCollection > &  trkColl,
const edm::ValueMap< float > &  trackTimes,
const edm::ValueMap< float > &  trackTimeResos 
) const

Definition at line 95 of file TransientTrackBuilder.cc.

References MillePedeFileConverter_cfg::e, mps_fire::i, edm::isNotFinite(), and hcalRecHitTable_cff::time.

97  {
98  vector<TransientTrack> ttVect;
99  ttVect.reserve((*trkColl).size());
100  for (unsigned int i = 0; i < (*trkColl).size(); i++) {
101  TrackRef ref(trkColl, i);
102  double time = trackTimes[ref];
103  double timeReso = trackTimeResos[ref];
104  timeReso = (timeReso > 1e-6 ? timeReso
105  : defaultInvalidTrackTimeReso); // make the error much larger than the BS time width
106  if (edm::isNotFinite(time)) {
107  time = 0.0;
108  timeReso = defaultInvalidTrackTimeReso;
109  }
110  ttVect.push_back(TransientTrack(ref, time, timeReso, theField, theTrackingGeometry));
111  }
112  return ttVect;
113 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
static constexpr float defaultInvalidTrackTimeReso

◆ build() [15/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::GsfTrackCollection > &  trkColl,
const edm::ValueMap< float > &  trackTimes,
const edm::ValueMap< float > &  trackTimeResos 
) const

Definition at line 115 of file TransientTrackBuilder.cc.

References MillePedeFileConverter_cfg::e, mps_fire::i, edm::isNotFinite(), and hcalRecHitTable_cff::time.

117  {
118  vector<TransientTrack> ttVect;
119  ttVect.reserve((*trkColl).size());
120  for (unsigned int i = 0; i < (*trkColl).size(); i++) {
121  GsfTrackRef ref(trkColl, i);
122  double time = trackTimes[ref];
123  double timeReso = trackTimeResos[ref];
124  timeReso = (timeReso > 1e-6 ? timeReso
125  : defaultInvalidTrackTimeReso); // make the error much larger than the BS time width
126  if (edm::isNotFinite(time)) {
127  time = 0.0;
128  timeReso = defaultInvalidTrackTimeReso;
129  }
130  ttVect.push_back(TransientTrack(new GsfTransientTrack(ref, time, timeReso, theField, theTrackingGeometry)));
131  }
132  return ttVect;
133 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
static constexpr float defaultInvalidTrackTimeReso

◆ build() [16/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< edm::View< reco::Track > > &  trkColl,
const edm::ValueMap< float > &  trackTimes,
const edm::ValueMap< float > &  trackTimeResos 
) const

Definition at line 135 of file TransientTrackBuilder.cc.

References MillePedeFileConverter_cfg::e, mps_fire::i, edm::isNotFinite(), and hcalRecHitTable_cff::time.

137  {
138  vector<TransientTrack> ttVect;
139  ttVect.reserve((*trkColl).size());
140  for (unsigned int i = 0; i < (*trkColl).size(); i++) {
141  const Track* trk = &(*trkColl)[i];
142  const GsfTrack* gsfTrack = dynamic_cast<const GsfTrack*>(trk);
143  if (gsfTrack) {
144  GsfTrackRef ref = RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>();
145  double time = trackTimes[ref];
146  double timeReso = trackTimeResos[ref];
147  timeReso = (timeReso > 1e-6 ? timeReso
148  : defaultInvalidTrackTimeReso); // make the error much larger than the BS time width
149  if (edm::isNotFinite(time)) {
150  time = 0.0;
151  timeReso = defaultInvalidTrackTimeReso;
152  }
153  ttVect.push_back(TransientTrack(new GsfTransientTrack(
154  RefToBase<Track>(trkColl, i).castTo<GsfTrackRef>(), time, timeReso, theField, theTrackingGeometry)));
155  } else { // gsf
156  TrackRef ref = RefToBase<Track>(trkColl, i).castTo<TrackRef>();
157  double time = trackTimes[ref];
158  double timeReso = trackTimeResos[ref];
159  timeReso = (timeReso > 1e-6 ? timeReso
160  : defaultInvalidTrackTimeReso); // make the error much larger than the BS time width
161  if (edm::isNotFinite(time)) {
162  time = 0.0;
163  timeReso = defaultInvalidTrackTimeReso;
164  }
165  ttVect.push_back(TransientTrack(
166  RefToBase<Track>(trkColl, i).castTo<TrackRef>(), time, timeReso, theField, theTrackingGeometry));
167  }
168  }
169  return ttVect;
170 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
const MagneticField * theField
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
static constexpr float defaultInvalidTrackTimeReso

◆ build() [17/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::TrackCollection > &  trkColl,
const reco::BeamSpot beamSpot 
) const

Definition at line 172 of file TransientTrackBuilder.cc.

References pwdgSkimBPark_cfi::beamSpot, newFWLiteAna::build, and mps_fire::i.

173  {
174  vector<TransientTrack> ttVect = build(trkColl);
175  for (unsigned int i = 0; i < ttVect.size(); i++) {
176  ttVect[i].setBeamSpot(beamSpot);
177  }
178  return ttVect;
179 }
reco::TransientTrack build(const reco::Track *p) const

◆ build() [18/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::GsfTrackCollection > &  trkColl,
const reco::BeamSpot beamSpot 
) const

Definition at line 181 of file TransientTrackBuilder.cc.

References pwdgSkimBPark_cfi::beamSpot, newFWLiteAna::build, and mps_fire::i.

182  {
183  vector<TransientTrack> ttVect = build(trkColl);
184  for (unsigned int i = 0; i < ttVect.size(); i++) {
185  ttVect[i].setBeamSpot(beamSpot);
186  }
187  return ttVect;
188 }
reco::TransientTrack build(const reco::Track *p) const

◆ build() [19/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< edm::View< reco::Track > > &  trkColl,
const reco::BeamSpot beamSpot 
) const

Definition at line 190 of file TransientTrackBuilder.cc.

References pwdgSkimBPark_cfi::beamSpot, newFWLiteAna::build, and mps_fire::i.

191  {
192  vector<TransientTrack> ttVect = build(trkColl);
193  for (unsigned int i = 0; i < ttVect.size(); i++) {
194  ttVect[i].setBeamSpot(beamSpot);
195  }
196  return ttVect;
197 }
reco::TransientTrack build(const reco::Track *p) const

◆ build() [20/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::TrackCollection > &  trkColl,
const reco::BeamSpot beamSpot,
const edm::ValueMap< float > &  trackTimes,
const edm::ValueMap< float > &  trackTimeResos 
) const

Definition at line 199 of file TransientTrackBuilder.cc.

References pwdgSkimBPark_cfi::beamSpot, newFWLiteAna::build, and mps_fire::i.

202  {
203  vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
204  for (unsigned int i = 0; i < ttVect.size(); i++) {
205  ttVect[i].setBeamSpot(beamSpot);
206  }
207  return ttVect;
208 }
reco::TransientTrack build(const reco::Track *p) const

◆ build() [21/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< reco::GsfTrackCollection > &  trkColl,
const reco::BeamSpot beamSpot,
const edm::ValueMap< float > &  trackTimes,
const edm::ValueMap< float > &  trackTimeResos 
) const

Definition at line 210 of file TransientTrackBuilder.cc.

References pwdgSkimBPark_cfi::beamSpot, newFWLiteAna::build, and mps_fire::i.

213  {
214  vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
215  for (unsigned int i = 0; i < ttVect.size(); i++) {
216  ttVect[i].setBeamSpot(beamSpot);
217  }
218  return ttVect;
219 }
reco::TransientTrack build(const reco::Track *p) const

◆ build() [22/23]

vector< TransientTrack > TransientTrackBuilder::build ( const edm::Handle< edm::View< reco::Track > > &  trkColl,
const reco::BeamSpot beamSpot,
const edm::ValueMap< float > &  trackTimes,
const edm::ValueMap< float > &  trackTimeResos 
) const

Definition at line 221 of file TransientTrackBuilder.cc.

References pwdgSkimBPark_cfi::beamSpot, newFWLiteAna::build, and mps_fire::i.

224  {
225  vector<TransientTrack> ttVect = build(trkColl, trackTimes, trackTimeResos);
226  for (unsigned int i = 0; i < ttVect.size(); i++) {
227  ttVect[i].setBeamSpot(beamSpot);
228  }
229  return ttVect;
230 }
reco::TransientTrack build(const reco::Track *p) const

◆ build() [23/23]

TransientTrack TransientTrackBuilder::build ( const FreeTrajectoryState fts) const

◆ field()

const MagneticField* TransientTrackBuilder::field ( ) const
inline

Definition at line 70 of file TransientTrackBuilder.h.

References theField.

70 { return theField; }
const MagneticField * theField

◆ trackingGeometry()

const edm::ESHandle<GlobalTrackingGeometry> TransientTrackBuilder::trackingGeometry ( ) const
inline

Definition at line 71 of file TransientTrackBuilder.h.

References theTrackingGeometry.

71 { return theTrackingGeometry; }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry

Member Data Documentation

◆ defaultInvalidTrackTimeReso

constexpr float TransientTrackBuilder::defaultInvalidTrackTimeReso = 0.350f
static

◆ theField

const MagneticField* TransientTrackBuilder::theField
private

Definition at line 75 of file TransientTrackBuilder.h.

Referenced by field().

◆ theTrackingGeometry

edm::ESHandle<GlobalTrackingGeometry> TransientTrackBuilder::theTrackingGeometry
private

Definition at line 76 of file TransientTrackBuilder.h.

Referenced by trackingGeometry().