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 
52 
53  unsigned int verbosity_;
54 
57 
60 
62 
64  {
73 
74  double ti_tr_min;
75  double ti_tr_max;
76 
77  void load(const edm::ParameterSet &ps)
78  {
79  x_cut_apply = ps.getParameter<bool> ("x_cut_apply");
80  x_cut_mean = ps.getParameter<double>("x_cut_mean");
81  x_cut_value = ps.getParameter<double>("x_cut_value");
82 
83  y_cut_apply = ps.getParameter<bool> ("y_cut_apply");
84  y_cut_mean = ps.getParameter<double>("y_cut_mean");
85  y_cut_value = ps.getParameter<double>("y_cut_value");
86 
87  xi_cut_apply = ps.getParameter<bool> ("xi_cut_apply");
88  xi_cut_mean = ps.getParameter<double>("xi_cut_mean");
89  xi_cut_value = ps.getParameter<double>("xi_cut_value");
90 
91  th_y_cut_apply = ps.getParameter<bool> ("th_y_cut_apply");
92  th_y_cut_mean = ps.getParameter<double>("th_y_cut_mean");
93  th_y_cut_value = ps.getParameter<double>("th_y_cut_value");
94 
95  ti_tr_min = ps.getParameter<double>("ti_tr_min");
96  ti_tr_max = ps.getParameter<double>("ti_tr_max");
97  }
98 
100  {
102 
103  desc.add<bool>("x_cut_apply", false)->setComment("whether to apply track-association cut in x");
104  desc.add<double>("x_cut_mean", 0E-6)->setComment("mean of track-association cut in x, mm");
105  desc.add<double>("x_cut_value", 800E-6)->setComment("threshold of track-association cut in x, mm");
106 
107  desc.add<bool>("y_cut_apply", false)->setComment("whether to apply track-association cut in y");
108  desc.add<double>("y_cut_mean", 0E-6)->setComment("mean of track-association cut in y, mm");
109  desc.add<double>("y_cut_value", 600E-6)->setComment("threshold of track-association cut in y, mm");
110 
111  desc.add<bool>("xi_cut_apply", true)->setComment("whether to apply track-association cut in xi");
112  desc.add<double>("xi_cut_mean", 0.)->setComment("mean of track-association cut in xi");
113  desc.add<double>("xi_cut_value", 0.013)->setComment("threshold of track-association cut in xi");
114 
115  desc.add<bool>("th_y_cut_apply", true)->setComment("whether to apply track-association cut in th_y");
116  desc.add<double>("th_y_cut_mean", 0E-6)->setComment("mean of track-association cut in th_y, rad");
117  desc.add<double>("th_y_cut_value", 20E-6)->setComment("threshold of track-association cut in th_y, rad");
118 
119  desc.add<double>("ti_tr_min", -1.)->setComment("minimum value for timing-tracking association cut");
120  desc.add<double>("ti_tr_max", +1.)->setComment("maximum value for timing-tracking association cut");
121 
122  return desc;
123  }
124  };
125 
126  std::map<unsigned int, AssociationCuts> association_cuts_; // map: arm -> AssociationCuts
127 
128  unsigned int max_n_timing_tracks_;
130 
132 
135 };
136 
137 //----------------------------------------------------------------------------------------------------
138 
140  tracksToken_ (consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagLocalTrackLite"))),
141 
142  pixelDiscardBXShiftedTracks_(iConfig.getParameter<bool>("pixelDiscardBXShiftedTracks")),
143 
144  lhcInfoLabel_ (iConfig.getParameter<std::string>("lhcInfoLabel")),
145  opticsLabel_ (iConfig.getParameter<std::string>("opticsLabel")),
146  verbosity_ (iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)),
147  doSingleRPReconstruction_ (iConfig.getParameter<bool>("doSingleRPReconstruction")),
148  doMultiRPReconstruction_ (iConfig.getParameter<bool>("doMultiRPReconstruction")),
149  singleRPReconstructionLabel_(iConfig.getParameter<std::string>("singleRPReconstructionLabel")),
150  multiRPReconstructionLabel_ (iConfig.getParameter<std::string>("multiRPReconstructionLabel")),
151 
152  localAngleXMin_ (iConfig.getParameter<double>("localAngleXMin")),
153  localAngleXMax_ (iConfig.getParameter<double>("localAngleXMax")),
154  localAngleYMin_ (iConfig.getParameter<double>("localAngleYMin")),
155  localAngleYMax_ (iConfig.getParameter<double>("localAngleYMax")),
156 
157  max_n_timing_tracks_ (iConfig.getParameter<unsigned int>("max_n_timing_tracks")),
158  default_time_ (iConfig.getParameter<double>("default_time")),
159 
160  algorithm_ (iConfig.getParameter<bool>("fitVtxY"), iConfig.getParameter<bool>("useImprovedInitialEstimate"), verbosity_),
162 {
163  for (const std::string &sector : { "45", "56" })
164  {
165  const unsigned int arm = (sector == "45") ? 0 : 1;
166  association_cuts_[arm].load(iConfig.getParameterSet("association_cuts_" + sector));
167  }
168 
170  produces<reco::ForwardProtonCollection>(singleRPReconstructionLabel_);
171 
173  produces<reco::ForwardProtonCollection>(multiRPReconstructionLabel_);
174 }
175 
176 //----------------------------------------------------------------------------------------------------
177 
179 {
181 
182  desc.add<edm::InputTag>("tagLocalTrackLite", edm::InputTag("ctppsLocalTrackLiteProducer"))
183  ->setComment("specification of the input lite-track collection");
184 
185  desc.add<bool>("pixelDiscardBXShiftedTracks", false)
186  ->setComment("whether to discard pixel tracks built from BX-shifted planes");
187 
188  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
189  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics record");
190 
191  desc.addUntracked<unsigned int>("verbosity", 0)->setComment("verbosity level");
192 
193  desc.add<bool>("doSingleRPReconstruction", true)
194  ->setComment("flag whether to apply single-RP reconstruction strategy");
195 
196  desc.add<bool>("doMultiRPReconstruction", true)
197  ->setComment("flag whether to apply multi-RP reconstruction strategy");
198 
199  desc.add<std::string>("singleRPReconstructionLabel", "singleRP")
200  ->setComment("output label for single-RP reconstruction products");
201 
202  desc.add<std::string>("multiRPReconstructionLabel", "multiRP")
203  ->setComment("output label for multi-RP reconstruction products");
204 
205  desc.add<double>("localAngleXMin", -0.03)->setComment("minimal accepted value of local horizontal angle (rad)");
206  desc.add<double>("localAngleXMax", +0.03)->setComment("maximal accepted value of local horizontal angle (rad)");
207  desc.add<double>("localAngleYMin", -0.04)->setComment("minimal accepted value of local vertical angle (rad)");
208  desc.add<double>("localAngleYMax", +0.04)->setComment("maximal accepted value of local vertical angle (rad)");
209 
210  for (const std::string &sector : { "45", "56" })
211  {
212  desc.add<edm::ParameterSetDescription>("association_cuts_" + sector, AssociationCuts::getDefaultParameters())
213  ->setComment("track-association cuts for sector " + sector);
214  }
215 
216  std::vector<edm::ParameterSet> config;
217 
218  desc.add<unsigned int>("max_n_timing_tracks", 5)->setComment("maximum number of timing tracks per RP");
219 
220  desc.add<double>("default_time", 0.)->setComment("proton time to be used when no timing information available");
221 
222  desc.add<bool>("fitVtxY", true)
223  ->setComment("for multi-RP reconstruction, flag whether y* should be free fit parameter");
224 
225  desc.add<bool>("useImprovedInitialEstimate", true)
226  ->setComment("for multi-RP reconstruction, flag whether a quadratic estimate of the initial point should be used");
227 
228  descriptions.add("ctppsProtons", desc);
229 }
230 
231 //----------------------------------------------------------------------------------------------------
232 
234 {
235  // get input
237  iEvent.getByToken(tracksToken_, hTracks);
238 
239  // book output
240  std::unique_ptr<reco::ForwardProtonCollection> pOutSingleRP(new reco::ForwardProtonCollection);
241  std::unique_ptr<reco::ForwardProtonCollection> pOutMultiRP(new reco::ForwardProtonCollection);
242 
243  // continue only if there is something to process
244  // NB: this avoids loading (possibly non-existing) conditions in workflows without proton data
245  if (!hTracks->empty())
246  {
247  // get conditions
248  edm::ESHandle<LHCInfo> hLHCInfo;
249  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
250 
252  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(opticsLabel_, hOpticalFunctions);
253 
255  iSetup.get<VeryForwardRealGeometryRecord>().get(hGeometry);
256 
257  // re-initialise algorithm upon crossing-angle change
258  if (opticsWatcher_.check(iSetup))
259  {
260  if (hOpticalFunctions->empty()) {
261  edm::LogInfo("CTPPSProtonProducer") << "No optical functions available, reconstruction disabled.";
263  opticsValid_ = false;
264  }
265  else {
266  algorithm_.init(*hOpticalFunctions);
267  opticsValid_ = true;
268  }
269  }
270 
271  // do reconstruction only if optics is valid
272  if (opticsValid_)
273  {
274  // prepare log
275  std::ostringstream ssLog;
276  if (verbosity_)
277  ssLog << "* input tracks:" << std::endl;
278 
279  // select tracks with small local angles, split them by LHC sector and tracker/timing RPs
280  std::map<unsigned int, std::vector<unsigned int>> trackingSelection, timingSelection;
281 
282  for (unsigned int idx = 0; idx < hTracks->size(); ++idx)
283  {
284  const auto& tr = hTracks->at(idx);
285 
286  if (tr.getTx() < localAngleXMin_ || tr.getTx() > localAngleXMax_
287  || tr.getTy() < localAngleYMin_ || tr.getTy() > localAngleYMax_)
288  continue;
289 
291  {
292  if (tr.getPixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::allShiftedPlanes ||
293  tr.getPixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::mixedPlanes)
294  continue;
295  }
296 
297  const CTPPSDetId rpId(tr.getRPId());
298 
299  if (verbosity_)
300  ssLog << "\t"
301  << "[" << idx << "] "
302  << tr.getRPId() << " (" << (rpId.arm()*100 + rpId.station()*10 + rpId.rp()) << "): "
303  << "x=" << tr.getX() << " +- " << tr.getXUnc() << " mm, "
304  << "y=" << tr.getY() << " +- " << tr.getYUnc() << " mm" << std::endl;
305 
306  const bool trackerRP = (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
307 
308  if (trackerRP)
309  trackingSelection[rpId.arm()].push_back(idx);
310  else
311  timingSelection[rpId.arm()].push_back(idx);
312  }
313 
314  // process each arm
315  for (const auto &arm_it : trackingSelection)
316  {
317  const auto &indices = arm_it.second;
318 
319  const auto &ac = association_cuts_[arm_it.first];
320 
321  // do single-RP reco if needed
322  std::map<unsigned int, reco::ForwardProton> singleRPResultsIndexed;
323  if (doSingleRPReconstruction_ || ac.xi_cut_apply || ac.th_y_cut_apply)
324  {
325  for (const auto &idx : indices)
326  {
327  if (verbosity_)
328  ssLog << std::endl << "* reconstruction from track " << idx << std::endl;
329 
330  singleRPResultsIndexed[idx] = algorithm_.reconstructFromSingleRP(CTPPSLocalTrackLiteRef(hTracks, idx), *hLHCInfo, ssLog);
331  }
332  }
333 
334  // check that exactly two tracking RPs are involved
335  // - 1 is insufficient for multi-RP reconstruction
336  // - PPS did not use more than 2 tracking RPs per arm -> algorithms are tuned to this
337  std::set<unsigned int> rpIds;
338  for (const auto &idx : indices)
339  rpIds.insert(hTracks->at(idx).getRPId());
340 
341  // do multi-RP reco if chosen
342  if (doMultiRPReconstruction_ && rpIds.size() == 2)
343  {
344  // find matching track pairs from different tracking RPs, ordered: i=near, j=far RP
345  std::vector<std::pair<unsigned int, unsigned int>> idx_pairs;
346  std::map<unsigned int, unsigned int> idx_pair_multiplicity;
347  for (const auto &i : indices)
348  {
349  for (const auto &j : indices)
350  {
351  const auto &tr_i = hTracks->at(i);
352  const auto &tr_j = hTracks->at(j);
353 
354  const double z_i = hGeometry->getRPTranslation(tr_i.getRPId()).z();
355  const double z_j = hGeometry->getRPTranslation(tr_j.getRPId()).z();
356 
357  const auto &pr_i = singleRPResultsIndexed[i];
358  const auto &pr_j = singleRPResultsIndexed[j];
359 
360  if (tr_i.getRPId() == tr_j.getRPId())
361  continue;
362 
363  if (std::abs(z_i) >= std::abs(z_j))
364  continue;
365 
366  bool matching = true;
367 
368  if (ac.x_cut_apply && std::abs(tr_i.getX() - tr_j.getX() - ac.x_cut_mean) > ac.x_cut_value)
369  matching = false;
370  else if (ac.y_cut_apply && std::abs(tr_i.getY() - tr_j.getY() - ac.y_cut_mean) > ac.y_cut_value)
371  matching = false;
372  else if (ac.xi_cut_apply && std::abs(pr_i.xi() - pr_j.xi() - ac.xi_cut_mean) > ac.xi_cut_value)
373  matching = false;
374  else if (ac.th_y_cut_apply && std::abs(pr_i.thetaY() - pr_j.thetaY() - ac.th_y_cut_mean) > ac.th_y_cut_value)
375  matching = false;
376 
377  if (!matching)
378  continue;
379 
380  idx_pairs.emplace_back(i, j);
381  idx_pair_multiplicity[i]++;
382  idx_pair_multiplicity[j]++;
383  }
384  }
385 
386  // evaluate track multiplicity in each timing RP
387  std::map<unsigned int, unsigned int> timing_RP_track_multiplicity;
388  for (const auto &ti : timingSelection[arm_it.first])
389  {
390  const auto &tr = hTracks->at(ti);
391  timing_RP_track_multiplicity[tr.getRPId()]++;
392  }
393 
394  // associate tracking-RP pairs with timing-RP tracks
395  std::map<unsigned int, std::vector<unsigned int>> matched_timing_track_indices;
396  std::map<unsigned int, unsigned int> matched_timing_track_multiplicity;
397  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx)
398  {
399  const auto &i = idx_pairs[pr_idx].first;
400  const auto &j = idx_pairs[pr_idx].second;
401 
402  // skip non-unique associations
403  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
404  continue;
405 
406  const auto &tr_i = hTracks->at(i);
407  const auto &tr_j = hTracks->at(j);
408 
409  const double z_i = hGeometry->getRPTranslation(tr_i.getRPId()).z();
410  const double z_j = hGeometry->getRPTranslation(tr_j.getRPId()).z();
411 
412  for (const auto &ti : timingSelection[arm_it.first])
413  {
414  const auto &tr_ti = hTracks->at(ti);
415 
416  // skip if timing RP saturated (high track multiplicity)
417  if (timing_RP_track_multiplicity[tr_ti.getRPId()] > max_n_timing_tracks_)
418  continue;
419 
420  // interpolation from tracking RPs
421  const double z_ti = - hGeometry->getRPTranslation(tr_ti.getRPId()).z(); // the minus sign fixes a bug in the diamond geometry
422  const double f_i = (z_ti - z_j) / (z_i - z_j), f_j = (z_i - z_ti) / (z_i - z_j);
423  const double x_inter = f_i * tr_i.getX() + f_j * tr_j.getX();
424  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();
425 
426  const double de_x = tr_ti.getX() - x_inter;
427  const double de_x_unc = sqrt(tr_ti.getXUnc()*tr_ti.getXUnc() + x_inter_unc_sq);
428  const double r = (de_x_unc > 0.) ? de_x / de_x_unc : 1E100;
429 
430  const bool matching = (ac.ti_tr_min < r && r < ac.ti_tr_max);
431 
432  if (verbosity_)
433  ssLog << "ti=" << ti << ", i=" << i << ", j=" << j
434  << " | z_ti=" << z_ti << ", z_i=" << z_i << ", z_j=" << z_j
435  << " | x_ti=" << tr_ti.getX() << ", x_inter=" << x_inter << ", de_x=" << de_x << ", de_x_unc=" << de_x_unc
436  << ", matching=" << matching << std::endl;
437 
438  if (!matching)
439  continue;
440 
441  matched_timing_track_indices[pr_idx].push_back(ti);
442  matched_timing_track_multiplicity[ti]++;
443  }
444  }
445 
446  // process associated tracks
447  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx)
448  {
449  const auto &i = idx_pairs[pr_idx].first;
450  const auto &j = idx_pairs[pr_idx].second;
451 
452  // skip non-unique associations of tracking-RP tracks
453  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
454  continue;
455 
456  if (verbosity_)
457  ssLog << std::endl << "* reconstruction from tracking-RP tracks: " << i << ", " << j << " and timing-RP tracks: ";
458 
459  // buffer contributing tracks
460  CTPPSLocalTrackLiteRefVector sel_tracks;
461  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, i));
462  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, j));
463 
464  CTPPSLocalTrackLiteRefVector sel_track_for_kin_reco = sel_tracks;
465 
466  // process timing-RP data
467  double sw=0., swt=0.;
468  for (const auto &ti : matched_timing_track_indices[pr_idx])
469  {
470  // skip non-unique associations of timing-RP tracks
471  if (matched_timing_track_multiplicity[ti] > 1)
472  continue;
473 
474  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, ti));
475 
476  if (verbosity_)
477  ssLog << ti << ", ";
478 
479  const auto &tr = hTracks->at(ti);
480  const double t_unc = tr.getTimeUnc();
481  const double w = (t_unc > 0.) ? 1./t_unc/t_unc : 1.;
482  sw += w;
483  swt += w * tr.getTime();
484  }
485 
486  float time = default_time_, time_unc = 0.;
487  if (sw > 0.)
488  {
489  time = swt / sw;
490  time_unc = 1. / sqrt(sw);
491  }
492 
493  if (verbosity_)
494  ssLog << std::endl << " time = " << time << " +- " << time_unc << std::endl;
495 
496  // process tracking-RP data
497  reco::ForwardProton proton = algorithm_.reconstructFromMultiRP(sel_track_for_kin_reco, *hLHCInfo, ssLog);
498 
499  // save combined output
500  proton.setContributingLocalTracks(sel_tracks);
501  proton.setTime(time);
502  proton.setTimeError(time_unc);
503 
504  pOutMultiRP->emplace_back(proton);
505  }
506  }
507 
508  // save single-RP results (un-indexed)
509  for (const auto &p : singleRPResultsIndexed)
510  pOutSingleRP->emplace_back(std::move(p.second));
511  }
512 
513  // dump log
514  if (verbosity_)
515  edm::LogInfo("CTPPSProtonProducer") << ssLog.str();
516  }
517  }
518 
519  // save output
521  iEvent.put(std::move(pOutSingleRP), singleRPReconstructionLabel_);
522 
524  iEvent.put(std::move(pOutMultiRP), multiRPReconstructionLabel_);
525 }
526 
527 //----------------------------------------------------------------------------------------------------
528 
530 
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_