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 };
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  produces<edm::ValueMap<float>>(t0Name);
72  produces<edm::ValueMap<float>>(sigmat0Name);
73  produces<edm::ValueMap<float>>(t0safeName);
74  produces<edm::ValueMap<float>>(sigmat0safeName);
75  produces<edm::ValueMap<float>>(probPiName);
76  produces<edm::ValueMap<float>>(probKName);
77  produces<edm::ValueMap<float>>(probPName);
78 }
79 
80 // Configuration descriptions
83  desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"))->setComment("Input tracks collection");
84  desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"))
85  ->setComment("Input ValueMap for track time at beamline");
86  desc.add<edm::InputTag>("tmtdSrc", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"))
87  ->setComment("Input ValueMap for track time at MTD");
88  desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"))
89  ->setComment("Input ValueMap for track time uncertainty at beamline");
90  desc.add<edm::InputTag>("sigmatmtdSrc", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"))
91  ->setComment("Input ValueMap for track time uncertainty at MTD");
92  desc.add<edm::InputTag>("pathLengthSrc", edm::InputTag("trackExtenderWithMTD:generalTrackPathLength"))
93  ->setComment("Input ValueMap for track path lengh from beamline to MTD");
94  desc.add<edm::InputTag>("pSrc", edm::InputTag("trackExtenderWithMTD:generalTrackp"))
95  ->setComment("Input ValueMap for track momentum magnitude (normally from refit with MTD hits)");
96  desc.add<edm::InputTag>("vtxsSrc", edm::InputTag("unsortedOfflinePrimaryVertices4DnoPID"))
97  ->setComment("Input primary vertex collection");
98  desc.add<double>("vtxMaxSigmaT", 0.025)
99  ->setComment("Maximum primary vertex time uncertainty for use in particle id [ns]");
100  desc.add<double>("maxDz", 0.1)
101  ->setComment("Maximum distance in z for track-primary vertex association for particle id [cm]");
102  desc.add<double>("maxDtSignificance", 5.0)
103  ->setComment(
104  "Maximum distance in time (normalized by uncertainty) for track-primary vertex association for particle id");
105  desc.add<double>("minProbHeavy", 0.75)
106  ->setComment("Minimum probability for a particle to be a kaon or proton before reassigning the timestamp");
107 
108  descriptions.add("tofPIDProducer", desc);
109 }
110 
111 template <class H, class T>
113  const edm::Handle<H>& handle,
114  const std::vector<T>& vec,
115  const std::string& name) const {
116  auto out = std::make_unique<edm::ValueMap<T>>();
118  filler.insert(handle, vec.begin(), vec.end());
119  filler.fill();
120  iEvent.put(std::move(out), name);
121 }
122 
124  constexpr double m_k = 0.493677; //[GeV]
125  constexpr double m_p = 0.9382720813; //[GeV]
126  constexpr double c_cm_ns = CLHEP::c_light * CLHEP::ns / CLHEP::cm; //[cm/ns]
127  constexpr double c_inv = 1.0 / c_cm_ns;
128 
130  ev.getByToken(tracksToken_, tracksH);
131  const auto& tracks = *tracksH;
132 
134  ev.getByToken(t0Token_, t0H);
135  const auto& t0In = *t0H;
136 
138  ev.getByToken(tmtdToken_, tmtdH);
139  const auto& tmtdIn = *tmtdH;
140 
142  ev.getByToken(sigmat0Token_, sigmat0H);
143  const auto& sigmat0In = *sigmat0H;
144 
146  ev.getByToken(sigmatmtdToken_, sigmatmtdH);
147  const auto& sigmatmtdIn = *sigmatmtdH;
148 
150  ev.getByToken(pathLengthToken_, pathLengthH);
151  const auto& pathLengthIn = *pathLengthH;
152 
154  ev.getByToken(pToken_, pH);
155  const auto& pIn = *pH;
156 
158  ev.getByToken(vtxsToken_, vtxsH);
159  const auto& vtxs = *vtxsH;
160 
161  //output value maps (PID probabilities and recalculated time at beamline)
162  std::vector<float> t0OutRaw;
163  std::vector<float> sigmat0OutRaw;
164  std::vector<float> t0safeOutRaw;
165  std::vector<float> sigmat0safeOutRaw;
166  std::vector<float> probPiOutRaw;
167  std::vector<float> probKOutRaw;
168  std::vector<float> probPOutRaw;
169 
170  //Do work here
171  for (unsigned int itrack = 0; itrack < tracks.size(); ++itrack) {
172  const reco::Track& track = tracks[itrack];
173  const reco::TrackRef trackref(tracksH, itrack);
174  float t0 = t0In[trackref];
175  float t0safe = t0;
176  float sigmat0safe = sigmat0In[trackref];
177  float sigmatmtd = sigmatmtdIn[trackref];
178  float sigmat0 = sigmatmtd;
179 
180  float prob_pi = -1.;
181  float prob_k = -1.;
182  float prob_p = -1.;
183 
184  if (sigmat0 > 0.) {
185  double rsigmazsq = 1. / track.dzError() / track.dzError();
186  double rsigmat = 1. / sigmatmtd;
187 
188  //find associated vertex
189  int vtxidx = -1;
190  int vtxidxmindz = -1;
191  int vtxidxminchisq = -1;
192  double mindz = maxDz_;
193  double minchisq = std::numeric_limits<double>::max();
194  //first try based on association weights, but keep track of closest in z and z-t as well
195  for (unsigned int ivtx = 0; ivtx < vtxs.size(); ++ivtx) {
196  const reco::Vertex& vtx = vtxs[ivtx];
197  float w = vtx.trackWeight(trackref);
198  if (w > 0.5) {
199  vtxidx = ivtx;
200  break;
201  }
202  double dz = std::abs(track.dz(vtx.position()));
203  if (dz < mindz) {
204  mindz = dz;
205  vtxidxmindz = ivtx;
206  }
207  if (vtx.tError() > 0. && vtx.tError() < vtxMaxSigmaT_) {
208  double dt = std::abs(t0 - vtx.t());
209  double dtsig = dt * rsigmat;
210  double chisq = dz * dz * rsigmazsq + dtsig * dtsig;
211  if (dz < maxDz_ && dtsig < maxDtSignificance_ && chisq < minchisq) {
212  minchisq = chisq;
213  vtxidxminchisq = ivtx;
214  }
215  }
216  }
217 
218  //if no vertex found based on association weights, fall back to closest in z or z-t
219  if (vtxidx < 0) {
220  //if closest vertex in z does not have valid time information, just use it,
221  //otherwise use the closest vertex in z-t plane with timing info, with a fallback to the closest in z
222  if (vtxidxmindz >= 0 && !(vtxs[vtxidxmindz].tError() > 0. && vtxs[vtxidxmindz].tError() < vtxMaxSigmaT_)) {
223  vtxidx = vtxidxmindz;
224  } else if (vtxidxminchisq >= 0) {
225  vtxidx = vtxidxminchisq;
226  } else if (vtxidxmindz >= 0) {
227  vtxidx = vtxidxmindz;
228  }
229  }
230 
231  //testing mass hypotheses only possible if there is an associated vertex with time information
232  if (vtxidx >= 0 && vtxs[vtxidx].tError() > 0. && vtxs[vtxidx].tError() < vtxMaxSigmaT_) {
233  //compute chisq in z-t plane for nominal vertex and mass hypothesis (pion)
234  const reco::Vertex& vtxnom = vtxs[vtxidx];
235  double dznom = std::abs(track.dz(vtxnom.position()));
236  double dtnom = std::abs(t0 - vtxnom.t());
237  double dtsignom = dtnom * rsigmat;
238  double chisqnom = dznom * dznom * rsigmazsq + dtsignom * dtsignom;
239 
240  //recompute t0 for alternate mass hypotheses
241  double t0_best = t0;
242 
243  //reliable match, revert to raw mtd time uncertainty
244  if (dtsignom < maxDtSignificance_) {
245  sigmat0safe = sigmatmtd;
246  }
247 
248  double tmtd = tmtdIn[trackref];
249  double pathlength = pathLengthIn[trackref];
250  double magp = pIn[trackref];
251 
252  double gammasq_k = 1. + magp * magp / m_k / m_k;
253  double beta_k = std::sqrt(1. - 1. / gammasq_k);
254  double t0_k = tmtd - pathlength / beta_k * c_inv;
255 
256  double gammasq_p = 1. + magp * magp / m_p / m_p;
257  double beta_p = std::sqrt(1. - 1. / gammasq_p);
258  double t0_p = tmtd - pathlength / beta_p * c_inv;
259 
260  double chisqmin = chisqnom;
261 
262  double chisqmin_pi = chisqnom;
263  double chisqmin_k = std::numeric_limits<double>::max();
264  double chisqmin_p = std::numeric_limits<double>::max();
265  //loop through vertices and check for better matches
266  for (const reco::Vertex& vtx : vtxs) {
267  if (!(vtx.tError() > 0. && vtx.tError() < vtxMaxSigmaT_)) {
268  continue;
269  }
270 
271  double dz = std::abs(track.dz(vtx.position()));
272  if (dz >= maxDz_) {
273  continue;
274  }
275 
276  double chisqdz = dz * dz * rsigmazsq;
277 
278  double dt_k = std::abs(t0_k - vtx.t());
279  double dtsig_k = dt_k * rsigmat;
280  double chisq_k = chisqdz + dtsig_k * dtsig_k;
281 
282  if (dtsig_k < maxDtSignificance_ && chisq_k < chisqmin_k) {
283  chisqmin_k = chisq_k;
284  }
285 
286  double dt_p = std::abs(t0_p - vtx.t());
287  double dtsig_p = dt_p * rsigmat;
288  double chisq_p = chisqdz + dtsig_p * dtsig_p;
289 
290  if (dtsig_p < maxDtSignificance_ && chisq_p < chisqmin_p) {
291  chisqmin_p = chisq_p;
292  }
293 
294  if (dtsig_k < maxDtSignificance_ && chisq_k < chisqmin) {
295  chisqmin = chisq_k;
296  t0_best = t0_k;
297  t0safe = t0_k;
298  sigmat0safe = sigmatmtd;
299  }
300  if (dtsig_p < maxDtSignificance_ && chisq_p < chisqmin) {
301  chisqmin = chisq_p;
302  t0_best = t0_p;
303  t0safe = t0_p;
304  sigmat0safe = sigmatmtd;
305  }
306  }
307 
308  //compute PID probabilities
309  //*TODO* deal with heavier nucleons and/or BSM case here?
310  double rawprob_pi = exp(-0.5 * chisqmin_pi);
311  double rawprob_k = exp(-0.5 * chisqmin_k);
312  double rawprob_p = exp(-0.5 * chisqmin_p);
313 
314  double normprob = 1. / (rawprob_pi + rawprob_k + rawprob_p);
315 
316  prob_pi = rawprob_pi * normprob;
317  prob_k = rawprob_k * normprob;
318  prob_p = rawprob_p * normprob;
319 
320  double prob_heavy = 1. - prob_pi;
321 
322  if (prob_heavy > minProbHeavy_) {
323  t0 = t0_best;
324  }
325  }
326  }
327 
328  t0OutRaw.push_back(t0);
329  sigmat0OutRaw.push_back(sigmat0);
330  t0safeOutRaw.push_back(t0safe);
331  sigmat0safeOutRaw.push_back(sigmat0safe);
332  probPiOutRaw.push_back(prob_pi);
333  probKOutRaw.push_back(prob_k);
334  probPOutRaw.push_back(prob_p);
335  }
336 
337  fillValueMap(ev, tracksH, t0OutRaw, t0Name);
338  fillValueMap(ev, tracksH, sigmat0OutRaw, sigmat0Name);
339  fillValueMap(ev, tracksH, t0safeOutRaw, t0safeName);
340  fillValueMap(ev, tracksH, sigmat0safeOutRaw, sigmat0safeName);
341  fillValueMap(ev, tracksH, probPiOutRaw, probPiName);
342  fillValueMap(ev, tracksH, probKOutRaw, probKName);
343  fillValueMap(ev, tracksH, probPOutRaw, probPName);
344 }
345 
346 //define this as a plug-in
void setComment(std::string const &value)
float dt
Definition: AMPTWrapper.h:136
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
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:525
edm::EDGetTokenT< edm::ValueMap< float > > sigmat0Token_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
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:113
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:19
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:84
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:596
double dzError() const
error on dz
Definition: TrackBase.h:725
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:129
#define constexpr
double t() const
t coordinate
Definition: Vertex.h:121