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