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:
24 
25  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
26 
27  template <class H, class T>
28  void fillValueMap(edm::Event& iEvent,
29  const edm::Handle<H>& handle,
30  const std::vector<T>& vec,
31  const std::string& name) const;
32 
33  void produce(edm::Event& ev, const edm::EventSetup& es) final;
34 
35 private:
36  static constexpr char t0Name[] = "t0";
37  static constexpr char sigmat0Name[] = "sigmat0";
38  static constexpr char t0safeName[] = "t0safe";
39  static constexpr char sigmat0safeName[] = "sigmat0safe";
40  static constexpr char probPiName[] = "probPi";
41  static constexpr char probKName[] = "probK";
42  static constexpr char probPName[] = "probP";
43 
56  const double vtxMaxSigmaT_;
57  const double maxDz_;
58  const double maxDtSignificance_;
59  const double minProbHeavy_;
60  const double fixedT0Error_;
61  const double probPion_;
62  const double probKaon_;
63  const double probProton_;
64  const double minTrackTimeQuality_;
65  const bool MVASel_;
66  const bool vertexReassignment_;
67 };
68 
70  : tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracksSrc"))),
71  t0Token_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0Src"))),
72  tmtdToken_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tmtdSrc"))),
73  sigmat0Token_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0Src"))),
74  sigmatmtdToken_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatmtdSrc"))),
75  tofkToken_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tofkSrc"))),
76  tofpToken_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tofpSrc"))),
77  sigmatofpiToken_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatofpiSrc"))),
78  sigmatofkToken_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatofkSrc"))),
79  sigmatofpToken_(consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatofpSrc"))),
80  vtxsToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vtxsSrc"))),
81  trackMTDTimeQualityToken_(
82  consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("trackMTDTimeQualityVMapTag"))),
83  vtxMaxSigmaT_(iConfig.getParameter<double>("vtxMaxSigmaT")),
84  maxDz_(iConfig.getParameter<double>("maxDz")),
85  maxDtSignificance_(iConfig.getParameter<double>("maxDtSignificance")),
86  minProbHeavy_(iConfig.getParameter<double>("minProbHeavy")),
87  fixedT0Error_(iConfig.getParameter<double>("fixedT0Error")),
88  probPion_(iConfig.getParameter<double>("probPion")),
89  probKaon_(iConfig.getParameter<double>("probKaon")),
90  probProton_(iConfig.getParameter<double>("probProton")),
91  minTrackTimeQuality_(iConfig.getParameter<double>("minTrackTimeQuality")),
92  MVASel_(iConfig.getParameter<bool>("MVASel")),
93  vertexReassignment_(iConfig.getParameter<bool>("vertexReassignment")) {
94  produces<edm::ValueMap<float>>(t0Name);
95  produces<edm::ValueMap<float>>(sigmat0Name);
96  produces<edm::ValueMap<float>>(t0safeName);
97  produces<edm::ValueMap<float>>(sigmat0safeName);
98  produces<edm::ValueMap<float>>(probPiName);
99  produces<edm::ValueMap<float>>(probKName);
100  produces<edm::ValueMap<float>>(probPName);
101 }
102 
103 // Configuration descriptions
106  desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"))->setComment("Input tracks collection");
107  desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"))
108  ->setComment("Input ValueMap for track time at beamline");
109  desc.add<edm::InputTag>("tmtdSrc", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"))
110  ->setComment("Input ValueMap for track time at MTD");
111  desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"))
112  ->setComment("Input ValueMap for track time uncertainty at beamline");
113  desc.add<edm::InputTag>("sigmatmtdSrc", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"))
114  ->setComment("Input ValueMap for track time uncertainty at MTD");
115  desc.add<edm::InputTag>("tofkSrc", edm::InputTag("trackExtenderWithMTD:generalTrackTofK"))
116  ->setComment("Input ValueMap for track tof as kaon");
117  desc.add<edm::InputTag>("tofpSrc", edm::InputTag("trackExtenderWithMTD:generalTrackTofP"))
118  ->setComment("Input ValueMap for track tof as proton");
119  desc.add<edm::InputTag>("sigmatofpiSrc", edm::InputTag("trackExtenderWithMTD:generalTrackSigmaTofPi"))
120  ->setComment("Input ValueMap for track sigma(tof) as pion");
121  desc.add<edm::InputTag>("sigmatofkSrc", edm::InputTag("trackExtenderWithMTD:generalTrackSigmaTofK"))
122  ->setComment("Input ValueMap for track sigma(tof) as kaon");
123  desc.add<edm::InputTag>("sigmatofpSrc", edm::InputTag("trackExtenderWithMTD:generalTrackSigmaTofP"))
124  ->setComment("Input ValueMap for track sigma(tof) as proton");
125  desc.add<edm::InputTag>("vtxsSrc", edm::InputTag("unsortedOfflinePrimaryVertices4DwithPID"))
126  ->setComment("Input primary vertex collection");
127  desc.add<edm::InputTag>("trackMTDTimeQualityVMapTag", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"))
128  ->setComment("Track MVA quality value");
129  desc.add<double>("vtxMaxSigmaT", 0.025)
130  ->setComment("Maximum primary vertex time uncertainty for use in particle id [ns]");
131  desc.add<double>("maxDz", 0.1)
132  ->setComment("Maximum distance in z for track-primary vertex association for particle id [cm]");
133  desc.add<double>("maxDtSignificance", 5.0)
134  ->setComment(
135  "Maximum distance in time (normalized by uncertainty) for track-primary vertex association for particle id");
136  desc.add<double>("minProbHeavy", 0.75)
137  ->setComment("Minimum probability for a particle to be a kaon or proton before reassigning the timestamp");
138  desc.add<double>("fixedT0Error", 0.)->setComment("Use a fixed T0 uncertainty [ns]");
139  desc.add<double>("probPion", 1.)->setComment("A priori probability pions");
140  desc.add<double>("probKaon", 1.)->setComment("A priori probability kaons");
141  desc.add<double>("probProton", 1.)->setComment("A priori probability for protons");
142  desc.add<double>("minTrackTimeQuality", 0.8)->setComment("Minimum MVA Quality selection on tracks");
143  desc.add<bool>("MVASel", false)->setComment("Use MVA Quality selection");
144  desc.add<bool>("vertexReassignment", true)->setComment("Track-vertex reassignment");
145 
146  descriptions.add("tofPIDProducer", desc);
147 }
148 
149 template <class H, class T>
151  const edm::Handle<H>& handle,
152  const std::vector<T>& vec,
153  const std::string& name) const {
154  auto out = std::make_unique<edm::ValueMap<T>>();
156  filler.insert(handle, vec.begin(), vec.end());
157  filler.fill();
158  iEvent.put(std::move(out), name);
159 }
160 
163  ev.getByToken(tracksToken_, tracksH);
164  const auto& tracks = *tracksH;
165 
166  const auto& t0In = ev.get(t0Token_);
167 
168  const auto& tmtdIn = ev.get(tmtdToken_);
169 
170  const auto& sigmat0In = ev.get(sigmat0Token_);
171 
172  const auto& sigmatmtdIn = ev.get(sigmatmtdToken_);
173 
174  const auto& tofkIn = ev.get(tofkToken_);
175 
176  const auto& tofpIn = ev.get(tofpToken_);
177 
178  const auto& sigmatofpiIn = ev.get(sigmatofpiToken_);
179 
180  const auto& sigmatofkIn = ev.get(sigmatofkToken_);
181 
182  const auto& sigmatofpIn = ev.get(sigmatofpToken_);
183 
184  const auto& vtxs = ev.get(vtxsToken_);
185 
186  const auto& trackMVAQualIn = ev.get(trackMTDTimeQualityToken_);
187 
188  //output value maps (PID probabilities and recalculated time at beamline)
189  std::vector<float> t0OutRaw;
190  std::vector<float> sigmat0OutRaw;
191  std::vector<float> t0safeOutRaw;
192  std::vector<float> sigmat0safeOutRaw;
193  std::vector<float> probPiOutRaw;
194  std::vector<float> probKOutRaw;
195  std::vector<float> probPOutRaw;
196 
197  //Do work here
198  for (unsigned int itrack = 0; itrack < tracks.size(); ++itrack) {
199  const reco::Track& track = tracks[itrack];
200  const reco::TrackRef trackref(tracksH, itrack);
201  float t0 = t0In[trackref];
202  float t0safe = t0;
203  float sigmat0safe = sigmat0In[trackref];
204  float sigmatmtd = (sigmatmtdIn[trackref] > 0. && fixedT0Error_ > 0.) ? fixedT0Error_ : sigmatmtdIn[trackref];
205  float sigmat0 = sigmatmtd;
206  float sigmatofpi = sigmatofpiIn[trackref];
207  float sigmatofk = sigmatofkIn[trackref];
208  float sigmatofp = sigmatofpIn[trackref];
209 
210  float prob_pi = -1.;
211  float prob_k = -1.;
212  float prob_p = -1.;
213 
214  float trackMVAQual = trackMVAQualIn[trackref];
215 
216  if (sigmat0 > 0. && (!MVASel_ || (MVASel_ && trackMVAQual >= minTrackTimeQuality_))) {
217  double rsigmazsq = 1. / track.dzError() / track.dzError();
218  std::array<double, 3> rsigmat = {{1. / std::sqrt(sigmatmtd * sigmatmtd + sigmatofpi * sigmatofpi),
219  1. / std::sqrt(sigmatmtd * sigmatmtd + sigmatofk * sigmatofk),
220  1. / std::sqrt(sigmatmtd * sigmatmtd + sigmatofp * sigmatofp)}};
221 
222  //find associated vertex
223  int vtxidx = -1;
224  int vtxidxmindz = -1;
225  int vtxidxminchisq = -1;
226  double mindz = maxDz_;
227  double minchisq = std::numeric_limits<double>::max();
228  //first try based on association weights, but keep track of closest in z and z-t as well
229  for (unsigned int ivtx = 0; ivtx < vtxs.size(); ++ivtx) {
230  const reco::Vertex& vtx = vtxs[ivtx];
231  float w = vtx.trackWeight(trackref);
232  if (w > 0.5) {
233  vtxidx = ivtx;
234  break;
235  }
236  double dz = std::abs(track.dz(vtx.position()));
237  if (dz < mindz) {
238  mindz = dz;
239  vtxidxmindz = ivtx;
240  }
241  if (vtx.tError() > 0. && vtx.tError() < vtxMaxSigmaT_) {
242  double dt = std::abs(t0 - vtx.t());
243  double dtsig = dt * rsigmat[0]; //pion hp. uncertainty
244  double chisq = dz * dz * rsigmazsq + dtsig * dtsig;
245  if (dz < maxDz_ && dtsig < maxDtSignificance_ && chisq < minchisq) {
246  minchisq = chisq;
247  vtxidxminchisq = ivtx;
248  }
249  }
250  }
251 
252  //if no vertex found based on association weights, fall back to closest in z or z-t
253  if (vtxidx < 0) {
254  //if closest vertex in z does not have valid time information, just use it,
255  //otherwise use the closest vertex in z-t plane with timing info, with a fallback to the closest in z
256  if (vtxidxmindz >= 0 && !(vtxs[vtxidxmindz].tError() > 0. && vtxs[vtxidxmindz].tError() < vtxMaxSigmaT_)) {
257  vtxidx = vtxidxmindz;
258  } else if (vtxidxminchisq >= 0) {
259  vtxidx = vtxidxminchisq;
260  } else if (vtxidxmindz >= 0) {
261  vtxidx = vtxidxmindz;
262  }
263  }
264 
265  //testing mass hypotheses only possible if there is an associated vertex with time information
266  if (vtxidx >= 0 && vtxs[vtxidx].tError() > 0. && vtxs[vtxidx].tError() < vtxMaxSigmaT_) {
267  //compute chisq in z-t plane for nominal vertex and mass hypothesis (pion)
268  const reco::Vertex& vtxnom = vtxs[vtxidx];
269  double dznom = std::abs(track.dz(vtxnom.position()));
270  double dtnom = std::abs(t0 - vtxnom.t());
271  double dtsignom = dtnom * rsigmat[0];
272  double chisqnom = dznom * dznom * rsigmazsq + dtsignom * dtsignom;
273 
274  //recompute t0 for alternate mass hypotheses
275  double t0_best = t0;
276 
277  //reliable match, revert to raw mtd time uncertainty + tof uncertainty for pion hp
278  if (dtsignom < maxDtSignificance_) {
279  sigmat0safe = 1. / rsigmat[0];
280  }
281 
282  double tmtd = tmtdIn[trackref];
283  double t0_k = tmtd - tofkIn[trackref];
284  double t0_p = tmtd - tofpIn[trackref];
285 
286  double chisqmin = chisqnom;
287 
288  double chisqmin_pi = chisqnom;
289  double chisqmin_k = std::numeric_limits<double>::max();
290  double chisqmin_p = std::numeric_limits<double>::max();
291  //loop through vertices and check for better matches
292  for (unsigned int ivtx = 0; ivtx < vtxs.size(); ++ivtx) {
293  const reco::Vertex& vtx = vtxs[ivtx];
294  if (!vertexReassignment_) {
295  if (ivtx != (unsigned int)vtxidx)
296  continue;
297  }
298  if (!(vtx.tError() > 0. && vtx.tError() < vtxMaxSigmaT_)) {
299  continue;
300  }
301 
302  double dz = std::abs(track.dz(vtx.position()));
303  if (dz >= maxDz_) {
304  continue;
305  }
306 
307  double chisqdz = dz * dz * rsigmazsq;
308 
309  double dt_k = std::abs(t0_k - vtx.t());
310  double dtsig_k = dt_k * rsigmat[1];
311  double chisq_k = chisqdz + dtsig_k * dtsig_k;
312 
313  if (dtsig_k < maxDtSignificance_ && chisq_k < chisqmin_k) {
314  chisqmin_k = chisq_k;
315  }
316 
317  double dt_p = std::abs(t0_p - vtx.t());
318  double dtsig_p = dt_p * rsigmat[2];
319  double chisq_p = chisqdz + dtsig_p * dtsig_p;
320 
321  if (dtsig_p < maxDtSignificance_ && chisq_p < chisqmin_p) {
322  chisqmin_p = chisq_p;
323  }
324 
325  if (dtsig_k < maxDtSignificance_ && chisq_k < chisqmin) {
326  chisqmin = chisq_k;
327  t0_best = t0_k;
328  t0safe = t0_k;
329  sigmat0safe = 1. / rsigmat[1];
330  }
331  if (dtsig_p < maxDtSignificance_ && chisq_p < chisqmin) {
332  chisqmin = chisq_p;
333  t0_best = t0_p;
334  t0safe = t0_p;
335  sigmat0safe = 1. / rsigmat[2];
336  }
337  }
338 
339  //compute PID probabilities
340  //*TODO* deal with heavier nucleons and/or BSM case here?
341  double rawprob_pi = probPion_ * exp(-0.5 * chisqmin_pi);
342  double rawprob_k = probKaon_ * exp(-0.5 * chisqmin_k);
343  double rawprob_p = probProton_ * exp(-0.5 * chisqmin_p);
344 
345  double normprob = 1. / (rawprob_pi + rawprob_k + rawprob_p);
346 
347  prob_pi = rawprob_pi * normprob;
348  prob_k = rawprob_k * normprob;
349  prob_p = rawprob_p * normprob;
350 
351  double prob_heavy = 1. - prob_pi;
352 
353  if (prob_heavy > minProbHeavy_) {
354  t0 = t0_best;
355  }
356  }
357  }
358 
359  t0OutRaw.push_back(t0);
360  sigmat0OutRaw.push_back(sigmat0);
361  t0safeOutRaw.push_back(t0safe);
362  sigmat0safeOutRaw.push_back(sigmat0safe);
363  probPiOutRaw.push_back(prob_pi);
364  probKOutRaw.push_back(prob_k);
365  probPOutRaw.push_back(prob_p);
366  }
367 
368  fillValueMap(ev, tracksH, t0OutRaw, t0Name);
369  fillValueMap(ev, tracksH, sigmat0OutRaw, sigmat0Name);
370  fillValueMap(ev, tracksH, t0safeOutRaw, t0safeName);
371  fillValueMap(ev, tracksH, sigmat0safeOutRaw, sigmat0safeName);
372  fillValueMap(ev, tracksH, probPiOutRaw, probPiName);
373  fillValueMap(ev, tracksH, probKOutRaw, probKName);
374  fillValueMap(ev, tracksH, probPOutRaw, probPName);
375 }
376 
377 //define this as a plug-in
float dt
Definition: AMPTWrapper.h:136
const bool MVASel_
edm::EDGetTokenT< reco::VertexCollection > vtxsToken_
T w() const
const Point & position() const
position
Definition: Vertex.h:128
edm::EDGetTokenT< edm::ValueMap< float > > sigmat0Token_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
const double minProbHeavy_
const double probPion_
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
const double fixedT0Error_
edm::EDGetTokenT< edm::ValueMap< float > > tofkToken_
static constexpr char t0Name[]
TOFPIDProducer(const ParameterSet &pset)
void produce(edm::Event &ev, const edm::EventSetup &es) final
static constexpr char t0safeName[]
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::ValueMap< float > > trackMTDTimeQualityToken_
T sqrt(T t)
Definition: SSEVec.h:19
const double vtxMaxSigmaT_
static constexpr char sigmat0safeName[]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const double probProton_
void fillValueMap(edm::Event &iEvent, const edm::Handle< H > &handle, const std::vector< T > &vec, const std::string &name) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< edm::ValueMap< float > > sigmatofpToken_
edm::EDGetTokenT< edm::ValueMap< float > > tmtdToken_
edm::EDGetTokenT< edm::ValueMap< float > > tofpToken_
edm::EDGetTokenT< edm::ValueMap< float > > sigmatofpiToken_
const bool vertexReassignment_
static constexpr char sigmat0Name[]
static constexpr char probKName[]
static constexpr char probPName[]
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const double minTrackTimeQuality_
double t() const
t coordinate
Definition: Vertex.h:136
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< edm::ValueMap< float > > t0Token_
fixed size matrix
HLT enums.
edm::EDGetTokenT< edm::ValueMap< float > > sigmatofkToken_
edm::EDGetTokenT< edm::ValueMap< float > > sigmatmtdToken_
const double probKaon_
const double maxDz_
const double maxDtSignificance_
static constexpr char probPiName[]
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
def move(src, dest)
Definition: eostools.py:511