CMS 3D CMS Logo

CTPPSProtonProducer.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  * Laurent Forthomme
5  ****************************************************************************/
6 
12 
14 
18 
21 
23 
26 
29 
32 
33 //----------------------------------------------------------------------------------------------------
34 
36 {
37  public:
38  explicit CTPPSProtonProducer(const edm::ParameterSet&);
39  ~CTPPSProtonProducer() override = default;
40 
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
42 
43  private:
44  void produce(edm::Event&, const edm::EventSetup&) override;
45 
47 
49 
50  unsigned int verbosity_;
51 
54 
57 
59 
61  {
63  double x_cut_value;
65  double y_cut_value;
67  double xi_cut_value;
70 
71  void load(const edm::ParameterSet &ps)
72  {
73  x_cut_apply = ps.getParameter<bool> ("x_cut_apply");
74  x_cut_value = ps.getParameter<double>("x_cut_value");
75  y_cut_apply = ps.getParameter<bool> ("y_cut_apply");
76  y_cut_value = ps.getParameter<double>("y_cut_value");
77  xi_cut_apply = ps.getParameter<bool> ("xi_cut_apply");
78  xi_cut_value = ps.getParameter<double>("xi_cut_value");
79  th_y_cut_apply = ps.getParameter<bool> ("th_y_cut_apply");
80  th_y_cut_value = ps.getParameter<double>("th_y_cut_value");
81  }
82 
84  {
86 
87  desc.add<bool>("x_cut_apply", false)->setComment("whether to apply track-association cut in x");
88  desc.add<double>("x_cut_value", 800E-6)->setComment("threshold of track-association cut in x, mm");
89  desc.add<bool>("y_cut_apply", false)->setComment("whether to apply track-association cut in y");
90  desc.add<double>("y_cut_value", 600E-6)->setComment("threshold of track-association cut in y, mm");
91  desc.add<bool>("xi_cut_apply", true)->setComment("whether to apply track-association cut in xi");
92  desc.add<double>("xi_cut_value", 0.013)->setComment("threshold of track-association cut in xi");
93  desc.add<bool>("th_y_cut_apply", true)->setComment("whether to apply track-association cut in th_y");
94  desc.add<double>("th_y_cut_value", 20E-6)->setComment("threshold of track-association cut in th_y, rad");
95 
96  return desc;
97  }
98  };
99 
100  std::map<unsigned int, AssociationCuts> association_cuts_; // map: arm -> AssociationCuts
101 
102  unsigned int max_n_timing_tracks_;
103 
105 
108 };
109 
110 //----------------------------------------------------------------------------------------------------
111 
113  tracksToken_ (consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagLocalTrackLite"))),
114  lhcInfoLabel_ (iConfig.getParameter<std::string>("lhcInfoLabel")),
115  verbosity_ (iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)),
116  doSingleRPReconstruction_ (iConfig.getParameter<bool>("doSingleRPReconstruction")),
117  doMultiRPReconstruction_ (iConfig.getParameter<bool>("doMultiRPReconstruction")),
118  singleRPReconstructionLabel_(iConfig.getParameter<std::string>("singleRPReconstructionLabel")),
119  multiRPReconstructionLabel_ (iConfig.getParameter<std::string>("multiRPReconstructionLabel")),
120 
121  localAngleXMin_ (iConfig.getParameter<double>("localAngleXMin")),
122  localAngleXMax_ (iConfig.getParameter<double>("localAngleXMax")),
123  localAngleYMin_ (iConfig.getParameter<double>("localAngleYMin")),
124  localAngleYMax_ (iConfig.getParameter<double>("localAngleYMax")),
125 
126  max_n_timing_tracks_ (iConfig.getParameter<unsigned int>("max_n_timing_tracks")),
127 
128  algorithm_ (iConfig.getParameter<bool>("fitVtxY"), iConfig.getParameter<bool>("useImprovedInitialEstimate"), verbosity_),
130 {
131  for (const std::string &sector : { "45", "56" })
132  {
133  const unsigned int arm = (sector == "45") ? 0 : 1;
134  association_cuts_[arm].load(iConfig.getParameterSet("association_cuts_" + sector));
135  }
136 
138  produces<reco::ForwardProtonCollection>(singleRPReconstructionLabel_);
139 
141  produces<reco::ForwardProtonCollection>(multiRPReconstructionLabel_);
142 }
143 
144 //----------------------------------------------------------------------------------------------------
145 
147 {
149 
150  desc.add<edm::InputTag>("tagLocalTrackLite", edm::InputTag("ctppsLocalTrackLiteProducer"))
151  ->setComment("specification of the input lite-track collection");
152 
153  desc.add<std::string>("lhcInfoLabel", "")
154  ->setComment("label of the LHCInfo record");
155 
156  desc.addUntracked<unsigned int>("verbosity", 0)->setComment("verbosity level");
157 
158  desc.add<bool>("doSingleRPReconstruction", true)
159  ->setComment("flag whether to apply single-RP reconstruction strategy");
160 
161  desc.add<bool>("doMultiRPReconstruction", true)
162  ->setComment("flag whether to apply multi-RP reconstruction strategy");
163 
164  desc.add<std::string>("singleRPReconstructionLabel", "singleRP")
165  ->setComment("output label for single-RP reconstruction products");
166 
167  desc.add<std::string>("multiRPReconstructionLabel", "multiRP")
168  ->setComment("output label for multi-RP reconstruction products");
169 
170  desc.add<double>("localAngleXMin", -0.03)->setComment("minimal accepted value of local horizontal angle (rad)");
171  desc.add<double>("localAngleXMax", +0.03)->setComment("maximal accepted value of local horizontal angle (rad)");
172  desc.add<double>("localAngleYMin", -0.04)->setComment("minimal accepted value of local vertical angle (rad)");
173  desc.add<double>("localAngleYMax", +0.04)->setComment("maximal accepted value of local vertical angle (rad)");
174 
175  for (const std::string &sector : { "45", "56" })
176  {
177  desc.add<edm::ParameterSetDescription>("association_cuts_" + sector, AssociationCuts::getDefaultParameters())
178  ->setComment("track-association cuts for sector " + sector);
179  }
180 
181  std::vector<edm::ParameterSet> config;
182 
183  desc.add<unsigned int>("max_n_timing_tracks", 5)->setComment("maximum number of timing tracks per RP");
184 
185  desc.add<bool>("fitVtxY", true)
186  ->setComment("for multi-RP reconstruction, flag whether y* should be free fit parameter");
187 
188  desc.add<bool>("useImprovedInitialEstimate", true)
189  ->setComment("for multi-RP reconstruction, flag whether a quadratic estimate of the initial point should be used");
190 
191  descriptions.add("ctppsProtons", desc);
192 }
193 
194 //----------------------------------------------------------------------------------------------------
195 
197 {
198  // get input
200  iEvent.getByToken(tracksToken_, hTracks);
201 
202  // book output
203  std::unique_ptr<reco::ForwardProtonCollection> pOutSingleRP(new reco::ForwardProtonCollection);
204  std::unique_ptr<reco::ForwardProtonCollection> pOutMultiRP(new reco::ForwardProtonCollection);
205 
206  // continue only if there is something to process
207  // NB: this avoids loading (possibly non-existing) conditions in workflows without proton data
208  if (!hTracks->empty())
209  {
210  // get conditions
211  edm::ESHandle<LHCInfo> hLHCInfo;
212  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
213 
215  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(hOpticalFunctions);
216 
218  iSetup.get<VeryForwardRealGeometryRecord>().get(hGeometry);
219 
220  // re-initialise algorithm upon crossing-angle change
221  if (opticsWatcher_.check(iSetup))
222  {
223  if (hOpticalFunctions->empty()) {
224  edm::LogInfo("CTPPSProtonProducer") << "No optical functions available, reconstruction disabled.";
226  opticsValid_ = false;
227  }
228  else {
229  algorithm_.init(*hOpticalFunctions);
230  opticsValid_ = true;
231  }
232  }
233 
234  // do reconstruction only if optics is valid
235  if (opticsValid_)
236  {
237  // prepare log
238  std::ostringstream ssLog;
239  if (verbosity_)
240  ssLog << "* input tracks:" << std::endl;
241 
242  // select tracks with small local angles, split them by LHC sector and tracker/timing RPs
243  std::map<unsigned int, std::vector<unsigned int>> trackingSelection, timingSelection;
244 
245  for (unsigned int idx = 0; idx < hTracks->size(); ++idx)
246  {
247  const auto& tr = hTracks->at(idx);
248 
249  if (tr.getTx() < localAngleXMin_ || tr.getTx() > localAngleXMax_
250  || tr.getTy() < localAngleYMin_ || tr.getTy() > localAngleYMax_)
251  continue;
252 
253  const CTPPSDetId rpId(tr.getRPId());
254 
255  if (verbosity_)
256  ssLog << "\t"
257  << "[" << idx << "] "
258  << tr.getRPId() << " (" << (rpId.arm()*100 + rpId.station()*10 + rpId.rp()) << "): "
259  << "x=" << tr.getX() << " +- " << tr.getXUnc() << " mm, "
260  << "y=" << tr.getY() << " +- " << tr.getYUnc() << " mm" << std::endl;
261 
262  const bool trackerRP = (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
263 
264  if (trackerRP)
265  trackingSelection[rpId.arm()].push_back(idx);
266  else
267  timingSelection[rpId.arm()].push_back(idx);
268  }
269 
270  // process each arm
271  for (const auto &arm_it : trackingSelection)
272  {
273  const auto &indices = arm_it.second;
274 
275  const auto &ac = association_cuts_[arm_it.first];
276 
277  // do single-RP reco if needed
278  std::map<unsigned int, reco::ForwardProton> singleRPResultsIndexed;
279  if (doSingleRPReconstruction_ || ac.xi_cut_apply || ac.th_y_cut_apply)
280  {
281  for (const auto &idx : indices)
282  {
283  if (verbosity_)
284  ssLog << std::endl << "* reconstruction from track " << idx << std::endl;
285 
286  singleRPResultsIndexed[idx] = algorithm_.reconstructFromSingleRP(CTPPSLocalTrackLiteRef(hTracks, idx), *hLHCInfo, ssLog);
287  }
288  }
289 
290  // check that exactly two tracking RPs are involved
291  // - 1 is insufficient for multi-RP reconstruction
292  // - PPS did not use more than 2 tracking RPs per arm -> algorithms are tuned to this
293  std::set<unsigned int> rpIds;
294  for (const auto &idx : indices)
295  rpIds.insert(hTracks->at(idx).getRPId());
296 
297  // do multi-RP reco if chosen
298  if (doMultiRPReconstruction_ && rpIds.size() == 2)
299  {
300  // find matching track pairs from different tracking RPs
301  std::vector<std::pair<unsigned int, unsigned int>> idx_pairs;
302  std::map<unsigned int, unsigned int> idx_pair_multiplicity;
303  for (const auto &i : indices)
304  {
305  for (const auto &j : indices)
306  {
307  if (j <= i)
308  continue;
309 
310  const auto &tr_i = hTracks->at(i);
311  const auto &tr_j = hTracks->at(j);
312 
313  const auto &pr_i = singleRPResultsIndexed[i];
314  const auto &pr_j = singleRPResultsIndexed[j];
315 
316  if (tr_i.getRPId() == tr_j.getRPId())
317  continue;
318 
319  bool matching = true;
320 
321  if (ac.x_cut_apply && std::abs(tr_i.getX() - tr_j.getX()) > ac.x_cut_value)
322  matching = false;
323  if (ac.y_cut_apply && std::abs(tr_i.getY() - tr_j.getY()) > ac.y_cut_value)
324  matching = false;
325  if (ac.xi_cut_apply && std::abs(pr_i.xi() - pr_j.xi()) > ac.xi_cut_value)
326  matching = false;
327  if (ac.th_y_cut_apply && std::abs(pr_i.thetaY() - pr_j.thetaY()) > ac.th_y_cut_value)
328  matching = false;
329 
330  if (!matching)
331  continue;
332 
333  idx_pairs.emplace_back(i, j);
334  idx_pair_multiplicity[i]++;
335  idx_pair_multiplicity[j]++;
336  }
337  }
338 
339  // evaluate track multiplicity in each timing RP
340  std::map<unsigned int, unsigned int> timing_RP_track_multiplicity;
341  for (const auto &ti : timingSelection[arm_it.first])
342  {
343  const auto &tr = hTracks->at(ti);
344  timing_RP_track_multiplicity[tr.getRPId()]++;
345  }
346 
347  // associate tracking-RP pairs with timing-RP tracks
348  std::map<unsigned int, std::vector<unsigned int>> matched_timing_track_indices;
349  std::map<unsigned int, unsigned int> matched_timing_track_multiplicity;
350  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx)
351  {
352  const auto &i = idx_pairs[pr_idx].first;
353  const auto &j = idx_pairs[pr_idx].second;
354 
355  // skip non-unique associations
356  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
357  continue;
358 
359  const auto &tr_i = hTracks->at(i);
360  const auto &tr_j = hTracks->at(j);
361 
362  const double z_i = hGeometry->getRPTranslation(tr_i.getRPId()).z();
363  const double z_j = hGeometry->getRPTranslation(tr_j.getRPId()).z();
364 
365  for (const auto &ti : timingSelection[arm_it.first])
366  {
367  const auto &tr_ti = hTracks->at(ti);
368 
369  // skip if timing RP saturated (high track multiplicity)
370  if (timing_RP_track_multiplicity[tr_ti.getRPId()] > max_n_timing_tracks_)
371  continue;
372 
373  // interpolation from tracking RPs
374  const double z_ti = - hGeometry->getRPTranslation(tr_ti.getRPId()).z(); // the minus sign fixes a bug in the diamond geometry
375  const double f_i = (z_ti - z_j) / (z_i - z_j), f_j = (z_i - z_ti) / (z_i - z_j);
376  const double x_inter = f_i * tr_i.getX() + f_j * tr_j.getX();
377  const double x_inter_unc_sq = f_i*f_i * tr_i.getXUnc()*tr_i.getXUnc() + f_j*f_j * tr_j.getXUnc()*tr_j.getXUnc();
378 
379  const double de_x = tr_ti.getX() - x_inter;
380  const double de_x_unc = sqrt(tr_ti.getXUnc()*tr_ti.getXUnc() + x_inter_unc_sq);
381 
382  const bool matching = (std::abs(de_x) <= de_x_unc);
383 
384  if (verbosity_)
385  ssLog << "ti=" << ti << ", i=" << i << ", j=" << j
386  << " | z_ti=" << z_ti << ", z_i=" << z_i << ", z_j=" << z_j
387  << " | x_ti=" << tr_ti.getX() << ", x_inter=" << x_inter << ", de_x=" << de_x << ", de_x_unc=" << de_x_unc
388  << ", matching=" << matching << std::endl;
389 
390  if (!matching)
391  continue;
392 
393  matched_timing_track_indices[pr_idx].push_back(ti);
394  matched_timing_track_multiplicity[ti]++;
395  }
396  }
397 
398  // process associated tracks
399  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx)
400  {
401  const auto &i = idx_pairs[pr_idx].first;
402  const auto &j = idx_pairs[pr_idx].second;
403 
404  // skip non-unique associations of tracking-RP tracks
405  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
406  continue;
407 
408  if (verbosity_)
409  ssLog << std::endl << "* reconstruction from tracking-RP tracks: " << i << ", " << j << " and timing-RP tracks: ";
410 
411  // process tracking-RP data
412  CTPPSLocalTrackLiteRefVector sel_tracks;
413  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, i));
414  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, j));
415  reco::ForwardProton proton = algorithm_.reconstructFromMultiRP(sel_tracks, *hLHCInfo, ssLog);
416 
417  // process timing-RP data
418  double sw=0., swt=0.;
419  for (const auto &ti : matched_timing_track_indices[pr_idx])
420  {
421  // skip non-unique associations of timing-RP tracks
422  if (matched_timing_track_multiplicity[ti] > 1)
423  continue;
424 
425  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, ti));
426 
427  if (verbosity_)
428  ssLog << ti << ", ";
429 
430  const auto &tr = hTracks->at(ti);
431  const double t_unc = tr.getTimeUnc();
432  const double w = (t_unc > 0.) ? 1./t_unc/t_unc : 1.;
433  sw += w;
434  swt += w * tr.getTime();
435  }
436 
437  float time = 0., time_unc = 0.;
438  if (sw > 0.)
439  {
440  time = swt / sw;
441  time_unc = 1. / sqrt(sw);
442  }
443 
444  if (verbosity_)
445  ssLog << std::endl << " time = " << time << " +- " << time_unc << std::endl;
446 
447  // save combined output
448  proton.setContributingLocalTracks(sel_tracks);
449  proton.setTime(time);
450  proton.setTimeError(time_unc);
451 
452  pOutMultiRP->emplace_back(proton);
453  }
454  }
455 
456  // save single-RP results (un-indexed)
457  for (const auto &p : singleRPResultsIndexed)
458  pOutSingleRP->emplace_back(std::move(p.second));
459  }
460 
461  // dump log
462  if (verbosity_)
463  edm::LogInfo("CTPPSProtonProducer") << ssLog.str();
464  }
465  }
466 
467  // save output
469  iEvent.put(std::move(pOutSingleRP), singleRPReconstructionLabel_);
470 
472  iEvent.put(std::move(pOutMultiRP), multiRPReconstructionLabel_);
473 }
474 
475 //----------------------------------------------------------------------------------------------------
476 
478 
T getParameter(std::string const &) const
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const double w
Definition: UKUtility.cc:23
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
reco::ForwardProton reconstructFromSingleRP(const CTPPSLocalTrackLiteRef &track, const LHCInfo &lhcInfo, std::ostream &os) const
run proton reconstruction using single-RP strategy
void load(const edm::ParameterSet &ps)
~CTPPSProtonProducer() override=default
unsigned int max_n_timing_tracks_
Event setup record containing the real (actual) geometry information.
void init(const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions)
config
Definition: looper.py:291
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setTimeError(float time_err)
T sqrt(T t)
Definition: SSEVec.h:18
reco::ForwardProton reconstructFromMultiRP(const CTPPSLocalTrackLiteRefVector &tracks, const LHCInfo &lhcInfo, std::ostream &os) const
run proton reconstruction using multiple-RP strategy
CTPPSProtonProducer(const edm::ParameterSet &)
edm::ESWatcher< CTPPSInterpolatedOpticsRcd > opticsWatcher_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Ref< CTPPSLocalTrackLiteCollection > CTPPSLocalTrackLiteRef
Persistent reference to a CTPPSLocalTrackLite.
std::map< unsigned int, AssociationCuts > association_cuts_
void setTime(float time)
static edm::ParameterSetDescription getDefaultParameters()
ParameterDescriptionBase * add(U const &iLabel, T const &value)
CLHEP::Hep3Vector getRPTranslation(unsigned int id) const
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
ParameterSet const & getParameterSet(std::string const &) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void produce(edm::Event &, const edm::EventSetup &) override
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLT enums.
T get() const
Definition: EventSetup.h:71
ProtonReconstructionAlgorithm algorithm_
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
void setContributingLocalTracks(const CTPPSLocalTrackLiteRefVector &v)
store the list of RP tracks that contributed to this global track
std::string singleRPReconstructionLabel_
def move(src, dest)
Definition: eostools.py:511
std::string multiRPReconstructionLabel_
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tracksToken_