CMS 3D CMS Logo

PATLeptonTimeLifeInfoProducer.cc
Go to the documentation of this file.
1 
16 
25 
36 
37 #include <cstring>
38 
39 template <typename T>
41 public:
44 
45  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
46  void produce(edm::Event&, const edm::EventSetup&) override;
47 
48 private:
49  //--- private utility methods
50  const reco::Track* getTrack(const T&);
53  static bool fitVertex(const std::vector<reco::TransientTrack>& transTrk, TransientVertex& transVtx) {
54  if (transTrk.size() < 2)
55  return false;
56  KalmanVertexFitter kvf(true);
57  transVtx = kvf.vertex(transTrk);
58  return transVtx.hasRefittedTracks() && transVtx.refittedTracks().size() == transTrk.size();
59  }
60 
61  //--- configuration parameters
66  int pvChoice_;
67 
69  //--- value map for TrackTimeLifeInfo (to be stored into the event)
71 };
72 
73 template <typename T>
75  : leptonsToken_(consumes<std::vector<T>>(cfg.getParameter<edm::InputTag>("src"))),
76  pvToken_(consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("pvSource"))),
77  transTrackBuilderToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
78  selector_(cfg.getParameter<std::string>("selection")),
79  pvChoice_(cfg.getParameter<int>("pvChoice")) {
80  produces<TrackTimeLifeInfoMap>();
81 }
82 
83 template <typename T>
85  // Get leptons
87  evt.getByToken(leptonsToken_, leptons);
88 
89  // Get the vertices
91  evt.getByToken(pvToken_, vertices);
92 
93  // Get transient track builder
94  const TransientTrackBuilder& transTrackBuilder = es.getData(transTrackBuilderToken_);
95 
96  std::vector<TrackTimeLifeInfo> infos;
97  infos.reserve(leptons->size());
98 
99  for (const auto& lepton : *leptons) {
101 
102  // Do nothing for lepton not passing selection
103  if (!selector_(lepton)) {
104  infos.push_back(info);
105  continue;
106  }
107  size_t pv_idx = 0;
108  if (pvChoice_ == useClosestInDz && getTrack(lepton) != nullptr) {
109  float dz_min = 999;
110  size_t vtx_idx = 0;
111  for (const auto& vtx : *vertices) {
112  float dz_tmp = std::abs(getTrack(lepton)->dz(vtx.position()));
113  if (dz_tmp < dz_min) {
114  dz_min = dz_tmp;
115  pv_idx = vtx_idx;
116  }
117  vtx_idx++;
118  }
119  }
120  const reco::Vertex& pv = !vertices->empty() ? (*vertices)[pv_idx] : reco::Vertex();
121 
122  // Obtain IP vector and set related info into lepton
123  produceAndFillIPInfo(lepton, transTrackBuilder, pv, info);
124 
125  // Fit SV and set related info for taus or do nothing for other lepton types
126  produceAndFillSVInfo(lepton, transTrackBuilder, pv, info);
127  infos.push_back(info);
128  } // end of lepton loop
129 
130  // Build the valuemap
131  auto infoMap = std::make_unique<TrackTimeLifeInfoMap>();
133  filler.insert(leptons, infos.begin(), infos.end());
134  filler.fill();
135 
136  // Store output into the event
137  evt.put(std::move(infoMap));
138 }
139 
140 template <>
142  return electron.gsfTrack().isNonnull() ? electron.gsfTrack().get() : nullptr;
143 }
144 
145 template <>
147  return muon.innerTrack().isNonnull() ? muon.innerTrack().get() : nullptr;
148 }
149 
150 template <>
152  const reco::Track* track = nullptr;
153  if (tau.leadChargedHadrCand().isNonnull())
154  track = tau.leadChargedHadrCand()->bestTrack();
155  return track;
156 }
157 
158 template <typename T>
160  const TransientTrackBuilder& transTrackBuilder,
161  const reco::Vertex& pv,
163  const reco::Track* track = getTrack(lepton);
164  if (track != nullptr) {
165  // Extrapolate track to the point closest to PV
166  reco::TransientTrack transTrack = transTrackBuilder.build(track);
167  AnalyticalImpactPointExtrapolator extrapolator(transTrack.field());
168  TrajectoryStateOnSurface closestState =
169  extrapolator.extrapolate(transTrack.impactPointState(), RecoVertex::convertPos(pv.position()));
170  if (!closestState.isValid()) {
171  edm::LogError("PATLeptonTimeLifeInfoProducer")
172  << "closestState not valid! From:\n"
173  << "transTrack.impactPointState():\n"
174  << transTrack.impactPointState() << "RecoVertex::convertPos(pv.position()):\n"
175  << RecoVertex::convertPos(pv.position());
176  return;
177  }
178  GlobalPoint pca = closestState.globalPosition();
179  GlobalError pca_cov = closestState.cartesianError().position();
180  GlobalVector ip_vec = GlobalVector(pca.x() - pv.x(), pca.y() - pv.y(), pca.z() - pv.z());
181  GlobalError ip_cov = pca_cov + GlobalError(pv.covariance());
182  VertexDistance3D pca_dist;
183  Measurement1D ip_mes = pca_dist.distance(pv, VertexState(pca, pca_cov));
184  if (ip_vec.dot(GlobalVector(lepton.px(), lepton.py(), lepton.pz())) < 0)
185  ip_mes = Measurement1D(-1. * ip_mes.value(), ip_mes.error());
186 
187  // Store Track and PCA info
188  info.setTrack(track);
189  info.setBField_z(transTrackBuilder.field()->inInverseGeV(GlobalPoint(track->vx(), track->vy(), track->vz())).z());
190  info.setPCA(pca, pca_cov);
191  info.setIP(ip_vec, ip_cov);
192  info.setIPLength(ip_mes);
193  }
194 }
195 
196 template <typename T>
198  const TransientTrackBuilder& transTrackBuilder,
199  const reco::Vertex& pv,
201 
202 template <>
204  const TransientTrackBuilder& transTrackBuilder,
205  const reco::Vertex& pv,
207  // Fit SV with tracks of charged tau decay products
208  int fitOK = 0;
209  if (tau.signalChargedHadrCands().size() + tau.signalLostTracks().size() > 1) {
210  // Get tracks from tau signal charged candidates
211  std::vector<reco::TransientTrack> transTrks;
212  TransientVertex transVtx;
213  for (const auto& cand : tau.signalChargedHadrCands()) {
214  if (cand.isNull())
215  continue;
216  const reco::Track* track = cand->bestTrack();
217  if (track != nullptr)
218  transTrks.push_back(transTrackBuilder.build(track));
219  }
220  for (const auto& cand : tau.signalLostTracks()) {
221  if (cand.isNull())
222  continue;
223  const reco::Track* track = cand->bestTrack();
224  if (track != nullptr)
225  transTrks.push_back(transTrackBuilder.build(track));
226  }
227  // Fit SV with KalmanVertexFitter
228  fitOK = fitVertex(transTrks, transVtx) ? 1 : -1;
229  if (fitOK > 0) {
230  reco::Vertex sv = transVtx;
231  // Get flight-length
232  // Full PV->SV flight vector with its covariance
233  GlobalVector flight_vec = GlobalVector(sv.x() - pv.x(), sv.y() - pv.y(), sv.z() - pv.z());
234  GlobalError flight_cov = transVtx.positionError() + GlobalError(pv.covariance());
235  //MB: can be taken from tau itself (but with different fit of PV and SV) as follows:
236  //tau.flightLength().mag2());
237  //tau.flightLengthSig();
238  VertexDistance3D sv_dist;
239  Measurement1D flightLength_mes = sv_dist.signedDistance(pv, sv, GlobalVector(tau.px(), tau.py(), tau.pz()));
240 
241  // Store SV info
242  info.setSV(sv);
243  info.setFlightVector(flight_vec, flight_cov);
244  info.setFlightLength(flightLength_mes);
245  }
246  }
247 }
248 
249 template <typename T>
251  // pat{Electron,Muon,Tau}TimeLifeInfoProducer
253 
254  std::string lepCollName;
255  if (typeid(T) == typeid(pat::Electron))
256  lepCollName = "slimmedElectrons";
257  else if (typeid(T) == typeid(pat::Muon))
258  lepCollName = "slimmedMuons";
259  else if (typeid(T) == typeid(pat::Tau))
260  lepCollName = "slimmedTaus";
261  desc.add<edm::InputTag>("src", edm::InputTag(lepCollName));
262  desc.add<edm::InputTag>("pvSource", edm::InputTag("offlineSlimmedPrimaryVertices"));
263  desc.add<std::string>("selection", "")->setComment("Selection required to produce and store time-life information");
264  desc.add<int>("pvChoice", useFront)
265  ->setComment(
266  "Define PV to compute IP: 0: first PV, 1: PV with the smallest dz of the tau leading track (default: " +
267  std::to_string(useFront) + ")");
268 
269  descriptions.addWithDefaultLabel(desc);
270 }
271 
reco::Vertex::Point convertPos(const GlobalPoint &p)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
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
T z() const
Definition: PV3DBase.h:61
GlobalError positionError() const
void produce(edm::Event &, const edm::EventSetup &) override
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
PATLeptonTimeLifeInfoProducer(const edm::ParameterSet &)
void produceAndFillSVInfo(const T &, const TransientTrackBuilder &, const reco::Vertex &, TrackTimeLifeInfo &)
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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
Log< level::Error, false > LogError
edm::EDGetTokenT< reco::VertexCollection > pvToken_
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
edm::EDGetTokenT< std::vector< T > > leptonsToken_
PATLeptonTimeLifeInfoProducer< pat::Electron > PATElectronTimeLifeInfoProducer
static std::string to_string(const XMLCh *ch)
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
reco::TransientTrack build(const reco::Track *p) const
PATLeptonTimeLifeInfoProducer< pat::Muon > PATMuonTimeLifeInfoProducer
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
static bool fitVertex(const std::vector< reco::TransientTrack > &transTrk, TransientVertex &transVtx)
const MagneticField * field() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Analysis-level tau class.
Definition: Tau.h:53
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
const MagneticField * field() const
const StringCutObjectSelector< T > selector_
bool hasRefittedTracks() const
void produceAndFillIPInfo(const T &, const TransientTrackBuilder &, const reco::Vertex &, TrackTimeLifeInfo &)
Analysis-level electron class.
Definition: Electron.h:51
double value() const
Definition: Measurement1D.h:25
fixed size matrix
const reco::Track * getTrack(const T &)
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double error() const
Definition: Measurement1D.h:27
long double T
Analysis-level muon class.
Definition: Muon.h:51
Structure to hold time-life information.
PATLeptonTimeLifeInfoProducer< pat::Tau > PATTauTimeLifeInfoProducer
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > transTrackBuilderToken_
def move(src, dest)
Definition: eostools.py:511
TrajectoryStateOnSurface impactPointState() const
Global3DVector GlobalVector
Definition: GlobalVector.h:10
std::vector< reco::TransientTrack > const & refittedTracks() const
Produces lepton life-time information.