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