CMS 3D CMS Logo

TrackTimeValueMapProducer.cc
Go to the documentation of this file.
1 // This producer assigns vertex times (with a specified resolution) to tracks.
2 // The times are produced as valuemaps associated to tracks, so the track
3 // dataformat doesn't need to be modified.
4 
7 
10 
12 
16 
32 
33 #include <memory>
34 #include <random>
35 
36 #include "CLHEP/Units/SystemOfUnits.h"
40 
42 public:
45 
46  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
47 
48 private:
49  // inputs
56  // tracking particle associators by order of preference
57  const std::vector<edm::EDGetTokenT<reco::TrackToTrackingParticleAssociator>> associators_;
58  // eta bounds
60  // options
61  std::vector<std::unique_ptr<const ResolutionModel>> resolutions_;
62  // functions
63  float extractTrackVertexTime(const TrackingParticle &, const reco::TransientTrack &) const;
64 };
65 
67 
68 namespace {
69  constexpr float m_pion = 139.57061e-3;
70  const std::string resolution("Resolution");
71 
72  template <typename ParticleType, typename T>
73  void writeValueMap(edm::Event &iEvent,
75  const std::vector<T> &values,
76  const std::string &label) {
77  std::unique_ptr<edm::ValueMap<T>> valMap(new edm::ValueMap<T>());
78  typename edm::ValueMap<T>::Filler filler(*valMap);
79  filler.insert(handle, values.begin(), values.end());
80  filler.fill();
81  iEvent.put(std::move(valMap), label);
82  }
83 } // namespace
84 
86  : tracks_(consumes<edm::View<reco::Track>>(conf.getParameter<edm::InputTag>("trackSrc"))),
87  tracksName_(conf.getParameter<edm::InputTag>("trackSrc").label()),
88  trackingParticles_(consumes<TrackingParticleCollection>(conf.getParameter<edm::InputTag>("trackingParticleSrc"))),
89  trackingVertices_(consumes<TrackingVertexCollection>(conf.getParameter<edm::InputTag>("trackingVertexSrc"))),
90  pileupSummaryInfo_(
91  consumes<std::vector<PileupSummaryInfo>>(conf.getParameter<edm::InputTag>("pileupSummaryInfo"))),
92  builderToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
93  associators_(edm::vector_transform(
94  conf.getParameter<std::vector<edm::InputTag>>("associators"),
95  [this](const edm::InputTag &tag) { return this->consumes<reco::TrackToTrackingParticleAssociator>(tag); })),
96  etaMin_(conf.getParameter<double>("etaMin")),
97  etaMax_(conf.getParameter<double>("etaMax")),
98  ptMin_(conf.getParameter<double>("ptMin")),
99  pMin_(conf.getParameter<double>("pMin")),
100  etaMaxForPtThreshold_(conf.getParameter<double>("etaMaxForPtThreshold")) {
101  // setup resolution models
102  const std::vector<edm::ParameterSet> &resos = conf.getParameterSetVector("resolutionModels");
103  for (const auto &reso : resos) {
104  const std::string &name = reso.getParameter<std::string>("modelName");
105  resolutions_.emplace_back(ResolutionModelFactory::get()->create(name, reso));
106 
107  // times and time resolutions for general tracks
108  produces<edm::ValueMap<float>>(tracksName_ + name);
109  produces<edm::ValueMap<float>>(tracksName_ + name + resolution);
110  }
111 }
112 
114  // get sim track associators
115  std::vector<edm::Handle<reco::TrackToTrackingParticleAssociator>> associators;
116  for (const auto &token : associators_) {
117  associators.emplace_back();
118  auto &back = associators.back();
119  evt.getByToken(token, back);
120  }
121 
122  std::vector<float> generalTrackTimes;
123 
124  // get track collections
125  edm::Handle<edm::View<reco::Track>> TrackCollectionH;
126  evt.getByToken(tracks_, TrackCollectionH);
127  const edm::View<reco::Track> &TrackCollection = *TrackCollectionH;
128 
129  // get tracking particle collections
131  evt.getByToken(trackingParticles_, TPCollectionH);
132 
134  evt.getByToken(pileupSummaryInfo_, pileupSummaryH);
135 
136  // transient track builder
137  auto const &builder = es.getData(builderToken_);
138 
139  // associate the reco tracks / gsf Tracks
140  std::vector<reco::RecoToSimCollection> associatedTracks;
141  associatedTracks.reserve(associators.size());
142  for (const auto &associator : associators) {
143  associatedTracks.emplace_back(associator->associateRecoToSim(TrackCollectionH, TPCollectionH));
144  }
145 
146  double sumSimTime = 0.;
147  double sumSimTimeSq = 0.;
148  int nsim = 0;
149  for (const PileupSummaryInfo &puinfo : *pileupSummaryH) {
150  if (puinfo.getBunchCrossing() == 0) {
151  for (const float &time : puinfo.getPU_times()) {
152  double simtime = time;
153  sumSimTime += simtime;
154  sumSimTimeSq += simtime * simtime;
155  ++nsim;
156  }
157  break;
158  }
159  }
160 
161  double meanSimTime = sumSimTime / double(nsim);
162  double varSimTime = sumSimTimeSq / double(nsim) - meanSimTime * meanSimTime;
163  double rmsSimTime = std::sqrt(std::max(0.1 * 0.1, varSimTime));
164  std::normal_distribution<float> gausSimTime(meanSimTime, rmsSimTime);
165 
166  // get event-based seed for RNG
167  unsigned int runNum_uint = static_cast<unsigned int>(evt.id().run());
168  unsigned int lumiNum_uint = static_cast<unsigned int>(evt.id().luminosityBlock());
169  unsigned int evNum_uint = static_cast<unsigned int>(evt.id().event());
170  unsigned int tkChi2_uint = uint32_t(TrackCollection.empty() ? 0 : TrackCollection.refAt(0)->chi2() / 0.01);
171  std::uint32_t seed = tkChi2_uint + (lumiNum_uint << 10) + (runNum_uint << 20) + evNum_uint;
172  std::mt19937 rng(seed);
173 
174  for (unsigned itk = 0; itk < TrackCollection.size(); ++itk) {
175  const auto tkref = TrackCollection.refAt(itk);
176  reco::RecoToSimCollection::const_iterator track_tps = associatedTracks.back().end();
177 
178  for (const auto &association : associatedTracks) {
179  track_tps = association.find(tkref);
180  if (track_tps != association.end())
181  break;
182  }
183 
184  if (track_tps != associatedTracks.back().end() && track_tps->val.size() == 1) {
185  reco::TransientTrack tt = builder.build(*tkref);
186  float time = extractTrackVertexTime(*track_tps->val[0].first, tt);
187  generalTrackTimes.push_back(time);
188  } else {
189  float rndtime = gausSimTime(rng);
190  generalTrackTimes.push_back(rndtime);
191  if (track_tps != associatedTracks.back().end() && track_tps->val.size() > 1) {
192  LogDebug("TooManyTracks") << "track matched to " << track_tps->val.size() << " tracking particles!"
193  << std::endl;
194  }
195  }
196  }
197 
198  for (const auto &reso : resolutions_) {
199  const std::string &name = reso->name();
200  std::vector<float> times, resos;
201 
202  times.reserve(TrackCollection.size());
203  resos.reserve(TrackCollection.size());
204 
205  for (unsigned i = 0; i < TrackCollection.size(); ++i) {
206  const reco::Track &tk = TrackCollection[i];
207  const float absEta = std::abs(tk.eta());
208  bool inAcceptance = absEta < etaMax_ && absEta >= etaMin_ && tk.p() > pMin_ &&
209  (absEta > etaMaxForPtThreshold_ || tk.pt() > ptMin_);
210  if (inAcceptance) {
211  const float resolution = reso->getTimeResolution(tk);
212  std::normal_distribution<float> gausGeneralTime(generalTrackTimes[i], resolution);
213  times.push_back(gausGeneralTime(rng));
214  resos.push_back(resolution);
215  } else {
216  times.push_back(0.0f);
217  resos.push_back(-1.);
218  }
219  }
220 
221  writeValueMap(evt, TrackCollectionH, times, tracksName_ + name);
222  writeValueMap(evt, TrackCollectionH, resos, tracksName_ + name + resolution);
223  }
224 }
225 
227  const reco::TransientTrack &tt) const {
228  int pdgid = tp.pdgId();
229  const auto &tvertex = tp.parentVertex();
230  math::XYZTLorentzVectorD result = tvertex->position();
231 
232  // account for secondary vertices...
233  if (tvertex->nSourceTracks() && tvertex->sourceTracks()[0]->pdgId() == pdgid) {
234  auto pvertex = tvertex->sourceTracks()[0]->parentVertex();
235  result = pvertex->position();
236  while (pvertex->nSourceTracks() && pvertex->sourceTracks()[0]->pdgId() == pdgid) {
237  pvertex = pvertex->sourceTracks()[0]->parentVertex();
238  result = pvertex->position();
239  }
240  }
241 
242  float time = result.T() * CLHEP::second;
243  // correct for time of flight from track reference position
244  GlobalPoint result_pos(result.x(), result.y(), result.z());
245  const auto &tkstate = tt.trajectoryStateClosestToPoint(result_pos);
246  float tkphi = tkstate.momentum().phi();
247  float tkz = tkstate.position().z();
248  float dphi = reco::deltaPhi(tkphi, tt.track().phi());
249  float dz = tkz - tt.track().vz();
250 
251  float radius = 100. * tt.track().pt() / (0.3 * tt.field()->inTesla(GlobalPoint(0, 0, 0)).z());
252  float pathlengthrphi = tt.track().charge() * dphi * radius;
253 
254  float pathlength = std::sqrt(pathlengthrphi * pathlengthrphi + dz * dz);
255  float p = tt.track().p();
256 
257  float speed = std::sqrt(1. / (1. + m_pion / p)) * CLHEP::c_light / CLHEP::cm; // speed in cm/ns
258  float dt = pathlength / speed;
259 
260  return time - dt;
261 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
float dt
Definition: AMPTWrapper.h:136
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::EDGetTokenT< TrackingVertexCollection > trackingVertices_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
def create(alignables, pedeDump, additionalData, outputFile, config)
double p() const
momentum vector magnitude
Definition: TrackBase.h:632
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::vector< std::unique_ptr< const ResolutionModel > > resolutions_
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:39
U second(std::pair< T, U > const &p)
double pt() const
track transverse momentum
Definition: TrackBase.h:638
char const * label
int iEvent
Definition: GenABIO.cc:224
Definition: TTTypes.h:54
std::tuple< layerClusterToCaloParticle, caloParticleToLayerCluster > association
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > builderToken_
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EventID id() const
Definition: EventBase.h:63
TrackTimeValueMapProducer(const edm::ParameterSet &)
const edm::EDGetTokenT< TrackingParticleCollection > trackingParticles_
RunNumber_t run() const
Definition: EventID.h:38
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:653
std::vector< TrackingVertex > TrackingVertexCollection
const std::vector< edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > > associators_
fixed size matrix
HLT enums.
const edm::EDGetTokenT< edm::View< reco::Track > > tracks_
Monte Carlo truth information used for tracking validation.
std::vector< TrackingParticle > TrackingParticleCollection
#define get
float extractTrackVertexTime(const TrackingParticle &, const reco::TransientTrack &) const
def move(src, dest)
Definition: eostools.py:511
EventNumber_t event() const
Definition: EventID.h:40
const edm::EDGetTokenT< std::vector< PileupSummaryInfo > > pileupSummaryInfo_
#define LogDebug(id)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override