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 dataformat doesn't
3 // need to be modified.
4 
7 
10 
12 
16 
32 
33 
34 #include <memory>
35 
37 #include "CLHEP/Units/SystemOfUnits.h"
39 #include "CLHEP/Random/RandGauss.h"
42 
45 
47 public:
50 
51  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
52 
53 private:
54  // inputs
60  // tracking particle associators by order of preference
61  const std::vector<edm::EDGetTokenT<reco::TrackToTrackingParticleAssociator> > associators_;
62  // eta bounds
64  // options
65  std::vector<std::unique_ptr<const ResolutionModel> > resolutions_;
66  // functions
68 };
69 
71 
72 namespace {
73  constexpr float m_pion = 139.57061e-3;
74  const std::string resolution("Resolution");
75 
76  template<typename ParticleType, typename T>
77  void writeValueMap(edm::Event &iEvent,
79  const std::vector<T> & values,
80  const std::string & label) {
81  std::unique_ptr<edm::ValueMap<T> > valMap(new edm::ValueMap<T>());
82  typename edm::ValueMap<T>::Filler filler(*valMap);
83  filler.insert(handle, values.begin(), values.end());
84  filler.fill();
85  iEvent.put(std::move(valMap), label);
86  }
87 }
88 
90  tracks_(consumes<edm::View<reco::Track> >( conf.getParameter<edm::InputTag>("trackSrc") ) ),
91  tracksName_(conf.getParameter<edm::InputTag>("trackSrc").label()),
92  trackingParticles_(consumes<TrackingParticleCollection>( conf.getParameter<edm::InputTag>("trackingParticleSrc") ) ),
93  trackingVertices_(consumes<TrackingVertexCollection>( conf.getParameter<edm::InputTag>("trackingVertexSrc") ) ),
94  pileupSummaryInfo_(consumes<std::vector<PileupSummaryInfo> >( conf.getParameter<edm::InputTag>("pileupSummaryInfo") ) ),
95  associators_( edm::vector_transform( conf.getParameter<std::vector<edm::InputTag> >("associators"), [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 {
102  // setup resolution models
103  const std::vector<edm::ParameterSet>& resos = conf.getParameterSetVector("resolutionModels");
104  for( const auto& reso : resos ) {
105  const std::string& name = reso.getParameter<std::string>("modelName");
106  ResolutionModel* resomod = ResolutionModelFactory::get()->create(name,reso);
107  resolutions_.emplace_back( resomod );
108 
109  // times and time resolutions for general tracks
110  produces<edm::ValueMap<float> >(tracksName_+name);
111  produces<edm::ValueMap<float> >(tracksName_+name+resolution);
112  }
113  // get RNG engine
115  if (!rng.isAvailable()){
116  throw cms::Exception("Configuration")
117  << "TrackTimeValueMapProducer::TrackTimeValueMapProducer() - RandomNumberGeneratorService is not present in configuration file.\n"
118  << "Add the service in the configuration file or remove the modules that require it.";
119  }
120 
121 }
122 
124  // get RNG engine
126  auto rng_engine = &(rng->getEngine(sid));
127 
128  // get sim track associators
129  std::vector<edm::Handle<reco::TrackToTrackingParticleAssociator> > associators;
130  for( const auto& token : associators_ ) {
131  associators.emplace_back();
132  auto& back = associators.back();
133  evt.getByToken(token,back);
134  }
135 
136  std::vector<float> generalTrackTimes;
137 
138  //get track collections
139  edm::Handle<edm::View<reco::Track> > TrackCollectionH;
140  evt.getByToken(tracks_, TrackCollectionH);
141  const edm::View<reco::Track>& TrackCollection = *TrackCollectionH;
142 
143  //get tracking particle collections
145  evt.getByToken(trackingParticles_, TPCollectionH);
146 
148  evt.getByToken(pileupSummaryInfo_, pileupSummaryH);
149 
150  //transient track builder
152  es.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
153 
154  // associate the reco tracks / gsf Tracks
155  std::vector<reco::RecoToSimCollection> associatedTracks;
156  for( auto associator : associators ) {
157  associatedTracks.emplace_back(associator->associateRecoToSim(TrackCollectionH, TPCollectionH));
158  }
159 
160 
161  double sumSimTime = 0.;
162  double sumSimTimeSq = 0.;
163  int nsim = 0;
164  for (const PileupSummaryInfo &puinfo : *pileupSummaryH) {
165  if (puinfo.getBunchCrossing() == 0) {
166  for (const float &time : puinfo.getPU_times()) {
167  double simtime = time;
168  sumSimTime += simtime;
169  sumSimTimeSq += simtime*simtime;
170  ++nsim;
171  }
172  break;
173  }
174  }
175 
176  double meanSimTime = sumSimTime/double(nsim);
177  double varSimTime = sumSimTimeSq/double(nsim) - meanSimTime*meanSimTime;
178  double rmsSimTime = std::sqrt(std::max(0.1*0.1,varSimTime));
179 
180  for( unsigned itk = 0; itk < TrackCollection.size(); ++itk ) {
181  const auto tkref = TrackCollection.refAt(itk);
182  reco::RecoToSimCollection::const_iterator track_tps = associatedTracks.back().end();
183 
184  for( const auto& association : associatedTracks ) {
185  track_tps = association.find(tkref);
186  if( track_tps != association.end() ) break;
187  }
188 
189  if (track_tps != associatedTracks.back().end() && track_tps->val.size() == 1) {
190  reco::TransientTrack tt = theB->build(*tkref);
191  float time = extractTrackVertexTime(*track_tps->val[0].first,tt);
192  generalTrackTimes.push_back(time);
193  }
194  else {
195  float rndtime = CLHEP::RandGauss::shoot(rng_engine, meanSimTime, rmsSimTime);
196  generalTrackTimes.push_back(rndtime);
197  if (track_tps != associatedTracks.back().end() && track_tps->val.size() > 1) {
198  LogDebug("TooManyTracks") << "track matched to " << track_tps->val.size() << " tracking particles!" << std::endl;
199  }
200  }
201  }
202 
203  for( const auto& reso : resolutions_ ) {
204  const std::string& name = reso->name();
205  std::vector<float> times, resos;
206 
207  times.reserve(TrackCollection.size());
208  resos.reserve(TrackCollection.size());
209 
210  for( unsigned i = 0; i < TrackCollection.size(); ++i ) {
211  const reco::Track& tk = TrackCollection[i];
212  const float absEta = std::abs(tk.eta());
213  bool inAcceptance = absEta < etaMax_ && absEta >= etaMin_ && tk.p()>pMin_ && (absEta>etaMaxForPtThreshold_ || tk.pt()>ptMin_);
214  if (inAcceptance) {
215  const float resolution = reso->getTimeResolution(tk);
216  times.push_back( CLHEP::RandGauss::shoot(rng_engine, generalTrackTimes[i], resolution) );
217  resos.push_back( resolution );
218  }
219  else {
220  times.push_back(0.0f);
221  resos.push_back(-1.);
222  }
223  }
224 
225  writeValueMap( evt, TrackCollectionH, times, tracksName_+name );
226  writeValueMap( evt, TrackCollectionH, resos, tracksName_+name+resolution );
227  }
228 }
229 
232  int pdgid = tp.pdgId();
233  const auto& tvertex = tp.parentVertex();
234  math::XYZTLorentzVectorD result = tvertex->position();
235 
236  // account for secondary vertices...
237  if( tvertex->nSourceTracks() && tvertex->sourceTracks()[0]->pdgId()==pdgid ) {
238  auto pvertex = tvertex->sourceTracks()[0]->parentVertex();
239  result = pvertex->position();
240  while( pvertex->nSourceTracks() && pvertex->sourceTracks()[0]->pdgId()==pdgid ) {
241  pvertex = pvertex->sourceTracks()[0]->parentVertex();
242  result = pvertex->position();
243  }
244  }
245 
246  float time = result.T()*CLHEP::second;
247  //correct for time of flight from track reference position
248  GlobalPoint result_pos(result.x(),result.y(),result.z());
249  const auto &tkstate = tt.trajectoryStateClosestToPoint(result_pos);
250  float tkphi = tkstate.momentum().phi();
251  float tkz = tkstate.position().z();
252  float dphi = reco::deltaPhi(tkphi,tt.track().phi());
253  float dz = tkz - tt.track().vz();
254 
255  float radius = 100.*tt.track().pt()/(0.3*tt.field()->inTesla(GlobalPoint(0,0,0)).z());
256  float pathlengthrphi = tt.track().charge()*dphi*radius;
257 
258  float pathlength = std::sqrt(pathlengthrphi*pathlengthrphi + dz*dz);
259  float p = tt.track().p();
260 
261  float speed = std::sqrt(1./(1.+m_pion/p))*CLHEP::c_light/CLHEP::cm; //speed in cm/ns
262  float dt = pathlength/speed;
263 
264  return time-dt;
265 }
#define LogDebug(id)
double p() const
momentum vector magnitude
Definition: TrackBase.h:615
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
float dt
Definition: AMPTWrapper.h:126
const edm::EDGetTokenT< TrackingVertexCollection > trackingVertices_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
std::vector< TrackingParticle > TrackingParticleCollection
friend struct const_iterator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int pdgId() const
PDG ID.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
reco::TransientTrack build(const reco::Track *p) const
size_type size() const
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
const MagneticField * field() const
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
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
#define constexpr
std::vector< std::unique_ptr< const ResolutionModel > > resolutions_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
RefToBase< value_type > refAt(size_type i) const
U second(std::pair< T, U > const &p)
int iEvent
Definition: GenABIO.cc:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:621
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float extractTrackVertexTime(const TrackingParticle &, const reco::TransientTrack &) const
double f[11][100]
const TrackingVertexRef & parentVertex() const
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
TrackTimeValueMapProducer(const edm::ParameterSet &)
const edm::EDGetTokenT< TrackingParticleCollection > trackingParticles_
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:669
const Track & track() const
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
std::vector< TrackingVertex > TrackingVertexCollection
const T & get() const
Definition: EventSetup.h:55
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.
int charge() const
track electric charge
Definition: TrackBase.h:567
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55
const edm::EDGetTokenT< std::vector< PileupSummaryInfo > > pileupSummaryInfo_