CMS 3D CMS Logo

TOFPIDProducer.cc
Go to the documentation of this file.
8 
15 #include "CLHEP/Units/GlobalPhysicalConstants.h"
16 #include "CLHEP/Units/GlobalSystemOfUnits.h"
17 
18 using namespace std;
19 using namespace edm;
20 
22  public:
23 
25 
26  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
27 
28  template <class H, class T>
29  void fillValueMap(edm::Event& iEvent, const edm::Handle<H>& handle, const std::vector<T>& vec, const std::string& name) const;
30 
31  void produce(edm::Event& ev, const edm::EventSetup& es) final;
32 
33  private:
34  static constexpr char t0Name[] = "t0";
35  static constexpr char sigmat0Name[] = "sigmat0";
36  static constexpr char t0safeName[] = "t0safe";
37  static constexpr char sigmat0safeName[] = "sigmat0safe";
38  static constexpr char probPiName[] = "probPi";
39  static constexpr char probKName[] = "probK";
40  static constexpr char probPName[] = "probP";
41 
50  double vtxMaxSigmaT_;
51  double maxDz_;
53  double minProbHeavy_;
54 };
55 
56 
57 
59  tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracksSrc"))),
60  t0Token_(consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("t0Src"))),
61  tmtdToken_(consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("tmtdSrc"))),
62  sigmat0Token_(consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("sigmat0Src"))),
63  sigmatmtdToken_(consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("sigmatmtdSrc"))),
64  pathLengthToken_(consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("pathLengthSrc"))),
65  pToken_(consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("pSrc"))),
66  vtxsToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vtxsSrc"))),
67  vtxMaxSigmaT_(iConfig.getParameter<double>("vtxMaxSigmaT")),
68  maxDz_(iConfig.getParameter<double>("maxDz")),
69  maxDtSignificance_(iConfig.getParameter<double>("maxDtSignificance")),
70  minProbHeavy_(iConfig.getParameter<double>("minProbHeavy"))
71  {
72  produces<edm::ValueMap<float> >(t0Name);
73  produces<edm::ValueMap<float> >(sigmat0Name);
74  produces<edm::ValueMap<float> >(t0safeName);
75  produces<edm::ValueMap<float> >(sigmat0safeName);
76  produces<edm::ValueMap<float> >(probPiName);
77  produces<edm::ValueMap<float> >(probKName);
78  produces<edm::ValueMap<float> >(probPName);
79 }
80 
81 // Configuration descriptions
84  desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"))->
85  setComment("Input tracks collection");
86  desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"))->
87  setComment("Input ValueMap for track time at beamline");
88  desc.add<edm::InputTag>("tmtdSrc", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"))->
89  setComment("Input ValueMap for track time at MTD");
90  desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"))->
91  setComment("Input ValueMap for track time uncertainty at beamline");
92  desc.add<edm::InputTag>("sigmatmtdSrc", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"))->
93  setComment("Input ValueMap for track time uncertainty at MTD");
94  desc.add<edm::InputTag>("pathLengthSrc", edm::InputTag("trackExtenderWithMTD:generalTrackPathLength"))->
95  setComment("Input ValueMap for track path lengh from beamline to MTD");
96  desc.add<edm::InputTag>("pSrc", edm::InputTag("trackExtenderWithMTD:generalTrackp"))->
97  setComment("Input ValueMap for track momentum magnitude (normally from refit with MTD hits)");
98  desc.add<edm::InputTag>("vtxsSrc", edm::InputTag("unsortedOfflinePrimaryVertices4DnoPID"))->
99  setComment("Input primary vertex collection");
100  desc.add<double>("vtxMaxSigmaT", 0.025)->
101  setComment("Maximum primary vertex time uncertainty for use in particle id [ns]");
102  desc.add<double>("maxDz", 0.1)->
103  setComment("Maximum distance in z for track-primary vertex association for particle id [cm]");
104  desc.add<double>("maxDtSignificance", 5.0)->
105  setComment("Maximum distance in time (normalized by uncertainty) for track-primary vertex association for particle id");
106  desc.add<double>("minProbHeavy", 0.75)->
107  setComment("Minimum probability for a particle to be a kaon or proton before reassigning the timestamp");
108 
109  descriptions.add("tofPIDProducer", desc);
110 }
111 
112 template <class H, class T>
113 void TOFPIDProducer::fillValueMap(edm::Event& iEvent, const edm::Handle<H>& handle, const std::vector<T>& vec, const std::string& name) const {
114  auto out = std::make_unique<edm::ValueMap<T>>();
116  filler.insert(handle, vec.begin(), vec.end());
117  filler.fill();
118  iEvent.put(std::move(out),name);
119 }
120 
122  constexpr double m_k = 0.493677; //[GeV]
123  constexpr double m_p = 0.9382720813; //[GeV]
124  constexpr double c_cm_ns = CLHEP::c_light*CLHEP::ns/CLHEP::cm; //[cm/ns]
125  constexpr double c_inv = 1.0/c_cm_ns;
126 
128  ev.getByToken(tracksToken_,tracksH);
129  const auto& tracks = *tracksH;
130 
132  ev.getByToken(t0Token_, t0H);
133  const auto &t0In = *t0H;
134 
136  ev.getByToken(tmtdToken_, tmtdH);
137  const auto &tmtdIn = *tmtdH;
138 
140  ev.getByToken(sigmat0Token_, sigmat0H);
141  const auto &sigmat0In = *sigmat0H;
142 
144  ev.getByToken(sigmatmtdToken_, sigmatmtdH);
145  const auto &sigmatmtdIn = *sigmatmtdH;
146 
147  edm::Handle<edm::ValueMap<float> > pathLengthH;
148  ev.getByToken(pathLengthToken_, pathLengthH);
149  const auto &pathLengthIn = *pathLengthH;
150 
152  ev.getByToken(pToken_, pH);
153  const auto &pIn = *pH;
154 
156  ev.getByToken(vtxsToken_,vtxsH);
157  const auto& vtxs = *vtxsH;
158 
159  //output value maps (PID probabilities and recalculated time at beamline)
160  std::vector<float> t0OutRaw;
161  std::vector<float> sigmat0OutRaw;
162  std::vector<float> t0safeOutRaw;
163  std::vector<float> sigmat0safeOutRaw;
164  std::vector<float> probPiOutRaw;
165  std::vector<float> probKOutRaw;
166  std::vector<float> probPOutRaw;
167 
168  //Do work here
169  for (unsigned int itrack = 0; itrack<tracks.size(); ++itrack) {
170  const reco::Track &track = tracks[itrack];
171  const reco::TrackRef trackref(tracksH,itrack);
172  float t0 = t0In[trackref];
173  float t0safe = t0;
174  float sigmat0safe = sigmat0In[trackref];
175  float sigmatmtd = sigmatmtdIn[trackref];
176  float sigmat0 = sigmatmtd;
177 
178  float prob_pi = -1.;
179  float prob_k = -1.;
180  float prob_p = -1.;
181 
182  if (sigmat0>0.) {
183 
184  double rsigmazsq = 1./track.dzError()/track.dzError();
185  double rsigmat = 1./sigmatmtd;
186 
187  //find associated vertex
188  int vtxidx = -1;
189  int vtxidxmindz = -1;
190  int vtxidxminchisq = -1;
191  double mindz = maxDz_;
192  double minchisq = std::numeric_limits<double>::max();
193  //first try based on association weights, but keep track of closest in z and z-t as well
194  for (unsigned int ivtx = 0; ivtx<vtxs.size(); ++ivtx) {
195  const reco::Vertex &vtx = vtxs[ivtx];
196  float w = vtx.trackWeight(trackref);
197  if (w>0.5) {
198  vtxidx = ivtx;
199  break;
200  }
201  double dz = std::abs(track.dz(vtx.position()));
202  if (dz<mindz) {
203  mindz = dz;
204  vtxidxmindz = ivtx;
205  }
206  if (vtx.tError()>0. && vtx.tError()<vtxMaxSigmaT_) {
207  double dt = std::abs(t0-vtx.t());
208  double dtsig = dt*rsigmat;
209  double chisq = dz*dz*rsigmazsq + dtsig*dtsig;
210  if (dz<maxDz_ && dtsig<maxDtSignificance_ && chisq<minchisq) {
211  minchisq = chisq;
212  vtxidxminchisq = ivtx;
213  }
214  }
215  }
216 
217  //if no vertex found based on association weights, fall back to closest in z or z-t
218  if (vtxidx<0) {
219  //if closest vertex in z does not have valid time information, just use it,
220  //otherwise use the closest vertex in z-t plane with timing info, with a fallback to the closest in z
221  if (vtxidxmindz>=0 && !(vtxs[vtxidxmindz].tError()>0. && vtxs[vtxidxmindz].tError()<vtxMaxSigmaT_)) {
222  vtxidx = vtxidxmindz;
223  }
224  else if (vtxidxminchisq>=0) {
225  vtxidx = vtxidxminchisq;
226  }
227  else if (vtxidxmindz>=0) {
228  vtxidx = vtxidxmindz;
229  }
230  }
231 
232  //testing mass hypotheses only possible if there is an associated vertex with time information
233  if (vtxidx>=0 && vtxs[vtxidx].tError()>0. && vtxs[vtxidx].tError()<vtxMaxSigmaT_) {
234  //compute chisq in z-t plane for nominal vertex and mass hypothesis (pion)
235  const reco::Vertex &vtxnom = vtxs[vtxidx];
236  double dznom = std::abs(track.dz(vtxnom.position()));
237  double dtnom = std::abs(t0 - vtxnom.t());
238  double dtsignom = dtnom*rsigmat;
239  double chisqnom = dznom*dznom*rsigmazsq + dtsignom*dtsignom;
240 
241  //recompute t0 for alternate mass hypotheses
242  double t0_best = t0;
243 
244  //reliable match, revert to raw mtd time uncertainty
245  if (dtsignom < maxDtSignificance_) {
246  sigmat0safe = sigmatmtd;
247  }
248 
249  double tmtd = tmtdIn[trackref];
250  double pathlength = pathLengthIn[trackref];
251  double magp = pIn[trackref];
252 
253  double gammasq_k = 1. + magp*magp/m_k/m_k;
254  double beta_k = std::sqrt(1.-1./gammasq_k);
255  double t0_k = tmtd - pathlength/beta_k*c_inv;
256 
257  double gammasq_p = 1. + magp*magp/m_p/m_p;
258  double beta_p = std::sqrt(1.-1./gammasq_p);
259  double t0_p = tmtd - pathlength/beta_p*c_inv;
260 
261  double chisqmin = chisqnom;
262 
263  double chisqmin_pi = chisqnom;
264  double chisqmin_k = std::numeric_limits<double>::max();
265  double chisqmin_p = std::numeric_limits<double>::max();
266  //loop through vertices and check for better matches
267  for (const reco::Vertex &vtx : vtxs) {
268  if (!(vtx.tError()>0. && vtx.tError()<vtxMaxSigmaT_)) {
269  continue;
270  }
271 
272  double dz = std::abs(track.dz(vtx.position()));
273  if (dz>=maxDz_) {
274  continue;
275  }
276 
277  double chisqdz = dz*dz*rsigmazsq;
278 
279  double dt_k = std::abs(t0_k - vtx.t());
280  double dtsig_k = dt_k*rsigmat;
281  double chisq_k = chisqdz + dtsig_k*dtsig_k;
282 
283  if (dtsig_k < maxDtSignificance_ && chisq_k<chisqmin_k) {
284  chisqmin_k = chisq_k;
285  }
286 
287  double dt_p = std::abs(t0_p - vtx.t());
288  double dtsig_p = dt_p*rsigmat;
289  double chisq_p = chisqdz + dtsig_p*dtsig_p;
290 
291  if (dtsig_p < maxDtSignificance_ && chisq_p<chisqmin_p) {
292  chisqmin_p = chisq_p;
293  }
294 
295  if (dtsig_k < maxDtSignificance_ && chisq_k<chisqmin) {
296  chisqmin = chisq_k;
297  t0_best = t0_k;
298  t0safe = t0_k;
299  sigmat0safe = sigmatmtd;
300  }
301  if (dtsig_p < maxDtSignificance_ && chisq_p<chisqmin) {
302  chisqmin = chisq_p;
303  t0_best = t0_p;
304  t0safe = t0_p;
305  sigmat0safe = sigmatmtd;
306  }
307 
308  }
309 
310  //compute PID probabilities
311  //*TODO* deal with heavier nucleons and/or BSM case here?
312  double rawprob_pi = exp(-0.5*chisqmin_pi);
313  double rawprob_k = exp(-0.5*chisqmin_k);
314  double rawprob_p = exp(-0.5*chisqmin_p);
315 
316  double normprob = 1./(rawprob_pi + rawprob_k + rawprob_p);
317 
318  prob_pi = rawprob_pi*normprob;
319  prob_k = rawprob_k*normprob;
320  prob_p = rawprob_p*normprob;
321 
322  double prob_heavy = 1.-prob_pi;
323 
324  if (prob_heavy>minProbHeavy_) {
325  t0 = t0_best;
326  }
327 
328  }
329 
330  }
331 
332  t0OutRaw.push_back(t0);
333  sigmat0OutRaw.push_back(sigmat0);
334  t0safeOutRaw.push_back(t0safe);
335  sigmat0safeOutRaw.push_back(sigmat0safe);
336  probPiOutRaw.push_back(prob_pi);
337  probKOutRaw.push_back(prob_k);
338  probPOutRaw.push_back(prob_p);
339  }
340 
341  fillValueMap(ev, tracksH, t0OutRaw, t0Name);
342  fillValueMap(ev, tracksH, sigmat0OutRaw, sigmat0Name);
343  fillValueMap(ev, tracksH, t0safeOutRaw, t0safeName);
344  fillValueMap(ev, tracksH, sigmat0safeOutRaw, sigmat0safeName);
345  fillValueMap(ev, tracksH, probPiOutRaw, probPiName);
346  fillValueMap(ev, tracksH, probKOutRaw, probKName);
347  fillValueMap(ev, tracksH, probPOutRaw, probPName);
348 
349 }
350 
351 //define this as a plug-in
float dt
Definition: AMPTWrapper.h:126
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
edm::EDGetTokenT< reco::VertexCollection > vtxsToken_
const double w
Definition: UKUtility.cc:23
static char probKName[]
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< edm::ValueMap< float > > sigmat0Token_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
edm::EDGetTokenT< edm::ValueMap< float > > pathLengthToken_
const Point & position() const
position
Definition: Vertex.h:109
double maxDtSignificance_
TOFPIDProducer(const ParameterSet &pset)
static char t0Name[]
void produce(edm::Event &ev, const edm::EventSetup &es) final
static char probPName[]
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:81
static char probPiName[]
edm::EDGetTokenT< edm::ValueMap< float > > tmtdToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static char sigmat0safeName[]
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:648
double dzError() const
error on dz
Definition: TrackBase.h:865
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static char t0safeName[]
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< edm::ValueMap< float > > t0Token_
fixed size matrix
HLT enums.
edm::EDGetTokenT< edm::ValueMap< float > > sigmatmtdToken_
edm::EDGetTokenT< edm::ValueMap< float > > pToken_
static char sigmat0Name[]
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
void fillValueMap(edm::Event &iEvent, const edm::Handle< H > &handle, const std::vector< T > &vec, const std::string &name) const
def move(src, dest)
Definition: eostools.py:511
double tError() const
error on t
Definition: Vertex.h:125
#define constexpr
double t() const
t coordinate
Definition: Vertex.h:117