CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes
PATLeptonTimeLifeInfoProducer< T > Class Template Reference

Produces lepton life-time information. More...

Inheritance diagram for PATLeptonTimeLifeInfoProducer< T >:
edm::stream::EDProducer<>

Public Member Functions

 PATLeptonTimeLifeInfoProducer (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~PATLeptonTimeLifeInfoProducer () override
 
- 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 Types

enum  PVChoice { useFront = 0, useClosestInDz }
 
using TrackTimeLifeInfoMap = edm::ValueMap< TrackTimeLifeInfo >
 

Private Member Functions

const reco::TrackgetTrack (const T &)
 
template<>
const reco::TrackgetTrack (const pat::Electron &electron)
 
template<>
const reco::TrackgetTrack (const pat::Muon &muon)
 
template<>
const reco::TrackgetTrack (const pat::Tau &tau)
 
void produceAndFillIPInfo (const T &, const TransientTrackBuilder &, const reco::Vertex &, TrackTimeLifeInfo &)
 
void produceAndFillSVInfo (const T &, const TransientTrackBuilder &, const reco::Vertex &, TrackTimeLifeInfo &)
 
template<>
void produceAndFillSVInfo (const pat::Tau &tau, const TransientTrackBuilder &transTrackBuilder, const reco::Vertex &pv, TrackTimeLifeInfo &info)
 

Static Private Member Functions

static bool fitVertex (const std::vector< reco::TransientTrack > &transTrk, TransientVertex &transVtx)
 

Private Attributes

edm::EDGetTokenT< std::vector< T > > leptonsToken_
 
int pvChoice_
 
edm::EDGetTokenT< reco::VertexCollectionpvToken_
 
const StringCutObjectSelector< Tselector_
 
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordtransTrackBuilderToken_
 

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

template<typename T>
class PATLeptonTimeLifeInfoProducer< T >

Produces lepton life-time information.

Author
Michal Bluj, NCBJ, Warsaw

Definition at line 39 of file PATLeptonTimeLifeInfoProducer.cc.

Member Typedef Documentation

◆ TrackTimeLifeInfoMap

Definition at line 69 of file PATLeptonTimeLifeInfoProducer.cc.

Member Enumeration Documentation

◆ PVChoice

template<typename T >
enum PATLeptonTimeLifeInfoProducer::PVChoice
private

Constructor & Destructor Documentation

◆ PATLeptonTimeLifeInfoProducer()

template<typename T >
PATLeptonTimeLifeInfoProducer< T >::PATLeptonTimeLifeInfoProducer ( const edm::ParameterSet cfg)
explicit

Definition at line 73 of file PATLeptonTimeLifeInfoProducer.cc.

74  : leptonsToken_(consumes<std::vector<T>>(cfg.getParameter<edm::InputTag>("src"))),
75  pvToken_(consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("pvSource"))),
76  transTrackBuilderToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
77  selector_(cfg.getParameter<std::string>("selection")),
78  pvChoice_(cfg.getParameter<int>("pvChoice")) {
79  produces<TrackTimeLifeInfoMap>();
80 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::EDGetTokenT< reco::VertexCollection > pvToken_
edm::EDGetTokenT< std::vector< T > > leptonsToken_
const StringCutObjectSelector< T > selector_
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transTrackBuilderToken_

◆ ~PATLeptonTimeLifeInfoProducer()

template<typename T >
PATLeptonTimeLifeInfoProducer< T >::~PATLeptonTimeLifeInfoProducer ( )
inlineoverride

Definition at line 42 of file PATLeptonTimeLifeInfoProducer.cc.

42 {};

Member Function Documentation

◆ fillDescriptions()

template<typename T >
void PATLeptonTimeLifeInfoProducer< T >::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 242 of file PATLeptonTimeLifeInfoProducer.cc.

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

242  {
243  // pat{Electron,Muon,Tau}TimeLifeInfoProducer
245 
246  std::string lepCollName;
247  if (typeid(T) == typeid(pat::Electron))
248  lepCollName = "slimmedElectrons";
249  else if (typeid(T) == typeid(pat::Muon))
250  lepCollName = "slimmedMuons";
251  else if (typeid(T) == typeid(pat::Tau))
252  lepCollName = "slimmedTaus";
253  desc.add<edm::InputTag>("src", edm::InputTag(lepCollName));
254  desc.add<edm::InputTag>("pvSource", edm::InputTag("offlineSlimmedPrimaryVertices"));
255  desc.add<std::string>("selection", "")->setComment("Selection required to produce and store time-life information");
256  desc.add<int>("pvChoice", useFront)
257  ->setComment(
258  "Define PV to compute IP: 0: first PV, 1: PV with the smallest dz of the tau leading track (default: " +
259  std::to_string(useFront) + ")");
260 
261  descriptions.addWithDefaultLabel(desc);
262 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
static std::string to_string(const XMLCh *ch)
Analysis-level tau class.
Definition: Tau.h:53
Analysis-level electron class.
Definition: Electron.h:51
long double T
Analysis-level muon class.
Definition: Muon.h:51

◆ fitVertex()

template<typename T >
static bool PATLeptonTimeLifeInfoProducer< T >::fitVertex ( const std::vector< reco::TransientTrack > &  transTrk,
TransientVertex transVtx 
)
inlinestaticprivate

Definition at line 52 of file PATLeptonTimeLifeInfoProducer.cc.

References TransientVertex::hasRefittedTracks(), TransientVertex::refittedTracks(), and KalmanVertexFitter::vertex().

52  {
53  if (transTrk.size() < 2)
54  return false;
55  KalmanVertexFitter kvf(true);
56  transVtx = kvf.vertex(transTrk);
57  return transVtx.hasRefittedTracks() && transVtx.refittedTracks().size() == transTrk.size();
58  }
bool hasRefittedTracks() const
std::vector< reco::TransientTrack > const & refittedTracks() const

◆ getTrack() [1/4]

template<typename T >
const reco::Track* PATLeptonTimeLifeInfoProducer< T >::getTrack ( const T )
private

◆ getTrack() [2/4]

template<>
const reco::Track * PATLeptonTimeLifeInfoProducer< pat::Electron >::getTrack ( const pat::Electron electron)
private

Definition at line 140 of file PATLeptonTimeLifeInfoProducer.cc.

References HPSPFTauProducerPuppi_cfi::electron.

140  {
141  return electron.gsfTrack().isNonnull() ? electron.gsfTrack().get() : nullptr;
142 }

◆ getTrack() [3/4]

template<>
const reco::Track * PATLeptonTimeLifeInfoProducer< pat::Muon >::getTrack ( const pat::Muon muon)
private

Definition at line 145 of file PATLeptonTimeLifeInfoProducer.cc.

145  {
146  return muon.innerTrack().isNonnull() ? muon.innerTrack().get() : nullptr;
147 }

◆ getTrack() [4/4]

template<>
const reco::Track * PATLeptonTimeLifeInfoProducer< pat::Tau >::getTrack ( const pat::Tau tau)
private

Definition at line 150 of file PATLeptonTimeLifeInfoProducer.cc.

References HLT_2024v10_cff::track.

150  {
151  const reco::Track* track = nullptr;
152  if (tau.leadChargedHadrCand().isNonnull())
153  track = tau.leadChargedHadrCand()->bestTrack();
154  return track;
155 }

◆ produce()

template<typename T >
void PATLeptonTimeLifeInfoProducer< T >::produce ( edm::Event evt,
const edm::EventSetup es 
)
override

Definition at line 83 of file PATLeptonTimeLifeInfoProducer.cc.

References funct::abs(), PVValHelper::dz, trigObjTnPSource_cfi::filler, edm::Event::getByToken(), edm::EventSetup::getData(), getTrack(), info(), submitPVValidationJobs::infos, HLT_2024v10_cff::leptons, eostools::move(), edm::Event::put(), HltBtagValidation_cff::Vertex, AlignmentTracksFromVertexSelector_cfi::vertices, and L1BJetProducer_cff::vtx.

83  {
84  // Get leptons
87 
88  // Get the vertices
91 
92  // Get transient track builder
93  const TransientTrackBuilder& transTrackBuilder = es.getData(transTrackBuilderToken_);
94 
95  std::vector<TrackTimeLifeInfo> infos;
96  infos.reserve(leptons->size());
97 
98  for (const auto& lepton : *leptons) {
100 
101  // Do nothing for lepton not passing selection
102  if (!selector_(lepton)) {
103  infos.push_back(info);
104  continue;
105  }
106  size_t pv_idx = 0;
107  if (pvChoice_ == useClosestInDz && getTrack(lepton) != nullptr) {
108  float dz_min = 999;
109  size_t vtx_idx = 0;
110  for (const auto& vtx : *vertices) {
111  float dz_tmp = std::abs(getTrack(lepton)->dz(vtx.position()));
112  if (dz_tmp < dz_min) {
113  dz_min = dz_tmp;
114  pv_idx = vtx_idx;
115  }
116  vtx_idx++;
117  }
118  }
119  const reco::Vertex& pv = !vertices->empty() ? (*vertices)[pv_idx] : reco::Vertex();
120 
121  // Obtain IP vector and set related info into lepton
122  produceAndFillIPInfo(lepton, transTrackBuilder, pv, info);
123 
124  // Fit SV and set related info for taus or do nothing for other lepton types
125  produceAndFillSVInfo(lepton, transTrackBuilder, pv, info);
126  infos.push_back(info);
127  } // end of lepton loop
128 
129  // Build the valuemap
130  auto infoMap = std::make_unique<TrackTimeLifeInfoMap>();
132  filler.insert(leptons, infos.begin(), infos.end());
133  filler.fill();
134 
135  // Store output into the event
136  evt.put(std::move(infoMap));
137 }
static const TGPicture * info(bool iBackgroundIsBlack)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
void produceAndFillSVInfo(const T &, const TransientTrackBuilder &, const reco::Vertex &, TrackTimeLifeInfo &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
edm::EDGetTokenT< reco::VertexCollection > pvToken_
edm::EDGetTokenT< std::vector< T > > leptonsToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
helper::Filler< ValueMap< T > > Filler
Definition: ValueMap.h:162
const StringCutObjectSelector< T > selector_
void produceAndFillIPInfo(const T &, const TransientTrackBuilder &, const reco::Vertex &, TrackTimeLifeInfo &)
const reco::Track * getTrack(const T &)
Structure to hold time-life information.
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transTrackBuilderToken_
def move(src, dest)
Definition: eostools.py:511

◆ produceAndFillIPInfo()

template<typename T >
void PATLeptonTimeLifeInfoProducer< T >::produceAndFillIPInfo ( const T lepton,
const TransientTrackBuilder transTrackBuilder,
const reco::Vertex pv,
TrackTimeLifeInfo info 
)
private

Definition at line 158 of file PATLeptonTimeLifeInfoProducer.cc.

References TransientTrackBuilder::build(), RecoVertex::convertPos(), Vector3DBase< T, FrameTag >::dot(), Measurement1D::error(), TransientTrackBuilder::field(), reco::TransientTrack::field(), getTrack(), reco::TransientTrack::impactPointState(), info(), MagneticField::inInverseGeV(), HLT_2024v10_cff::track, Measurement1D::value(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

161  {
162  const reco::Track* track = getTrack(lepton);
163  if (track != nullptr) {
164  info.setTrack(track);
165  info.setBField_z(transTrackBuilder.field()->inInverseGeV(GlobalPoint(track->vx(), track->vy(), track->vz())).z());
166 
167  // Extrapolate track to the point closest to PV
168  reco::TransientTrack transTrack = transTrackBuilder.build(track);
169  AnalyticalImpactPointExtrapolator extrapolator(transTrack.field());
170  TrajectoryStateOnSurface closestState =
171  extrapolator.extrapolate(transTrack.impactPointState(), RecoVertex::convertPos(pv.position()));
172  GlobalPoint pca = closestState.globalPosition();
173  GlobalError pca_cov = closestState.cartesianError().position();
174  GlobalVector ip_vec = GlobalVector(pca.x() - pv.x(), pca.y() - pv.y(), pca.z() - pv.z());
175  GlobalError ip_cov = pca_cov + GlobalError(pv.covariance());
176  VertexDistance3D pca_dist;
177  Measurement1D ip_mes = pca_dist.distance(pv, VertexState(pca, pca_cov));
178  if (ip_vec.dot(GlobalVector(lepton.px(), lepton.py(), lepton.pz())) < 0)
179  ip_mes = Measurement1D(-1. * ip_mes.value(), ip_mes.error());
180 
181  // Store PCA info
182  info.setPCA(pca, pca_cov);
183  info.setIP(ip_vec, ip_cov);
184  info.setIPLength(ip_mes);
185  }
186 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
static const TGPicture * info(bool iBackgroundIsBlack)
T z() const
Definition: PV3DBase.h:61
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:36
GlobalErrorBase< double, ErrorMatrixTag > GlobalError
Definition: GlobalError.h:13
reco::TransientTrack build(const reco::Track *p) const
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const MagneticField * field() const
const MagneticField * field() const
double value() const
Definition: Measurement1D.h:25
const reco::Track * getTrack(const T &)
double error() const
Definition: Measurement1D.h:27
TrajectoryStateOnSurface impactPointState() const
Global3DVector GlobalVector
Definition: GlobalVector.h:10

◆ produceAndFillSVInfo() [1/2]

template<typename T >
void PATLeptonTimeLifeInfoProducer< T >::produceAndFillSVInfo ( const T lepton,
const TransientTrackBuilder transTrackBuilder,
const reco::Vertex pv,
TrackTimeLifeInfo info 
)
private

Definition at line 189 of file PATLeptonTimeLifeInfoProducer.cc.

192  {}

◆ produceAndFillSVInfo() [2/2]

template<>
void PATLeptonTimeLifeInfoProducer< pat::Tau >::produceAndFillSVInfo ( const pat::Tau tau,
const TransientTrackBuilder transTrackBuilder,
const reco::Vertex pv,
TrackTimeLifeInfo info 
)
private

Definition at line 195 of file PATLeptonTimeLifeInfoProducer.cc.

References TransientTrackBuilder::build(), info(), TransientVertex::positionError(), pfDeepBoostedJetPreprocessParams_cfi::sv, and HLT_2024v10_cff::track.

198  {
199  // Fit SV with tracks of charged tau decay products
200  int fitOK = 0;
201  if (tau.signalChargedHadrCands().size() + tau.signalLostTracks().size() > 1) {
202  // Get tracks from tau signal charged candidates
203  std::vector<reco::TransientTrack> transTrks;
204  TransientVertex transVtx;
205  for (const auto& cand : tau.signalChargedHadrCands()) {
206  if (cand.isNull())
207  continue;
208  const reco::Track* track = cand->bestTrack();
209  if (track != nullptr)
210  transTrks.push_back(transTrackBuilder.build(track));
211  }
212  for (const auto& cand : tau.signalLostTracks()) {
213  if (cand.isNull())
214  continue;
215  const reco::Track* track = cand->bestTrack();
216  if (track != nullptr)
217  transTrks.push_back(transTrackBuilder.build(track));
218  }
219  // Fit SV with KalmanVertexFitter
220  fitOK = fitVertex(transTrks, transVtx) ? 1 : -1;
221  if (fitOK > 0) {
222  reco::Vertex sv = transVtx;
223  // Get flight-length
224  // Full PV->SV flight vector with its covariance
225  GlobalVector flight_vec = GlobalVector(sv.x() - pv.x(), sv.y() - pv.y(), sv.z() - pv.z());
226  GlobalError flight_cov = transVtx.positionError() + GlobalError(pv.covariance());
227  //MB: can be taken from tau itself (but with different fit of PV and SV) as follows:
228  //tau.flightLength().mag2());
229  //tau.flightLengthSig();
230  VertexDistance3D sv_dist;
231  Measurement1D flightLength_mes = sv_dist.signedDistance(pv, sv, GlobalVector(tau.px(), tau.py(), tau.pz()));
232 
233  // Store SV info
234  info.setSV(sv);
235  info.setFlightVector(flight_vec, flight_cov);
236  info.setFlightLength(flightLength_mes);
237  }
238  }
239 }
static const TGPicture * info(bool iBackgroundIsBlack)
GlobalError positionError() const
GlobalErrorBase< double, ErrorMatrixTag > GlobalError
Definition: GlobalError.h:13
reco::TransientTrack build(const reco::Track *p) const
static bool fitVertex(const std::vector< reco::TransientTrack > &transTrk, TransientVertex &transVtx)
Global3DVector GlobalVector
Definition: GlobalVector.h:10

Member Data Documentation

◆ leptonsToken_

template<typename T >
edm::EDGetTokenT<std::vector<T> > PATLeptonTimeLifeInfoProducer< T >::leptonsToken_
private

Definition at line 61 of file PATLeptonTimeLifeInfoProducer.cc.

◆ pvChoice_

template<typename T >
int PATLeptonTimeLifeInfoProducer< T >::pvChoice_
private

Definition at line 65 of file PATLeptonTimeLifeInfoProducer.cc.

◆ pvToken_

template<typename T >
edm::EDGetTokenT<reco::VertexCollection> PATLeptonTimeLifeInfoProducer< T >::pvToken_
private

Definition at line 62 of file PATLeptonTimeLifeInfoProducer.cc.

◆ selector_

template<typename T >
const StringCutObjectSelector<T> PATLeptonTimeLifeInfoProducer< T >::selector_
private

Definition at line 64 of file PATLeptonTimeLifeInfoProducer.cc.

◆ transTrackBuilderToken_

template<typename T >
edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> PATLeptonTimeLifeInfoProducer< T >::transTrackBuilderToken_
private

Definition at line 63 of file PATLeptonTimeLifeInfoProducer.cc.