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