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 
14 
16 
20 
23 
25 
33 
36 
39 
42 
43 //----------------------------------------------------------------------------------------------------
44 
46 public:
47  explicit CTPPSProtonProducer(const edm::ParameterSet &);
48  ~CTPPSProtonProducer() override = default;
49 
50  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
51 
52 private:
53  void produce(edm::Event &, const edm::EventSetup &) override;
54 
56 
58 
64  const bool useNewLHCInfo_;
65 
66  unsigned int verbosity_;
67 
70 
73 
75 
76  unsigned int max_n_timing_tracks_;
77  double default_time_;
78 
80 
83 
90 };
91 
92 //----------------------------------------------------------------------------------------------------
93 
95  : tracksToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagLocalTrackLite"))),
96 
97  pixelDiscardBXShiftedTracks_(iConfig.getParameter<bool>("pixelDiscardBXShiftedTracks")),
98 
99  lhcInfoPerLSLabel_(iConfig.getParameter<std::string>("lhcInfoPerLSLabel")),
100  lhcInfoPerFillLabel_(iConfig.getParameter<std::string>("lhcInfoPerFillLabel")),
101  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
102  opticsLabel_(iConfig.getParameter<std::string>("opticsLabel")),
103  ppsAssociationCutsLabel_(iConfig.getParameter<std::string>("ppsAssociationCutsLabel")),
104  useNewLHCInfo_(iConfig.getParameter<bool>("useNewLHCInfo")),
105  verbosity_(iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)),
106  doSingleRPReconstruction_(iConfig.getParameter<bool>("doSingleRPReconstruction")),
107  doMultiRPReconstruction_(iConfig.getParameter<bool>("doMultiRPReconstruction")),
108  singleRPReconstructionLabel_(iConfig.getParameter<std::string>("singleRPReconstructionLabel")),
109  multiRPReconstructionLabel_(iConfig.getParameter<std::string>("multiRPReconstructionLabel")),
110 
111  localAngleXMin_(iConfig.getParameter<double>("localAngleXMin")),
112  localAngleXMax_(iConfig.getParameter<double>("localAngleXMax")),
113  localAngleYMin_(iConfig.getParameter<double>("localAngleYMin")),
114  localAngleYMax_(iConfig.getParameter<double>("localAngleYMax")),
115 
116  max_n_timing_tracks_(iConfig.getParameter<unsigned int>("max_n_timing_tracks")),
117  default_time_(iConfig.getParameter<double>("default_time")),
118 
119  algorithm_(iConfig.getParameter<bool>("fitVtxY"),
120  iConfig.getParameter<bool>("useImprovedInitialEstimate"),
121  iConfig.getParameter<std::string>("multiRPAlgorithm"),
122  verbosity_),
123  opticsValid_(false),
124  lhcInfoPerLSToken_(esConsumes<LHCInfoPerLS, LHCInfoPerLSRcd>(edm::ESInputTag("", lhcInfoPerLSLabel_))),
125  lhcInfoPerFillToken_(esConsumes<LHCInfoPerFill, LHCInfoPerFillRcd>(edm::ESInputTag("", lhcInfoPerFillLabel_))),
126  lhcInfoToken_(esConsumes<LHCInfo, LHCInfoRcd>(edm::ESInputTag("", lhcInfoLabel_))),
128  edm::ESInputTag("", opticsLabel_))),
130  ppsAssociationCutsToken_(
131  esConsumes<PPSAssociationCuts, PPSAssociationCutsRcd>(edm::ESInputTag("", ppsAssociationCutsLabel_))) {
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<bool>("pixelDiscardBXShiftedTracks", false)
148  ->setComment("whether to discard pixel tracks built from BX-shifted planes");
149 
150  desc.add<std::string>("lhcInfoPerFillLabel", "")->setComment("label of the LHCInfoPerFill record");
151  desc.add<std::string>("lhcInfoPerLSLabel", "")->setComment("label of the LHCInfoPerLS record");
152  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
153  desc.add<bool>("useNewLHCInfo", false)->setComment("flag whether to use new LHCInfoPer* records or old LHCInfo");
154 
155  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics record");
156  desc.add<std::string>("ppsAssociationCutsLabel", "")->setComment("label of the association cuts record");
157 
158  desc.addUntracked<unsigned int>("verbosity", 0)->setComment("verbosity level");
159 
160  desc.add<bool>("doSingleRPReconstruction", true)
161  ->setComment("flag whether to apply single-RP reconstruction strategy");
162 
163  desc.add<bool>("doMultiRPReconstruction", true)->setComment("flag whether to apply multi-RP reconstruction strategy");
164 
165  desc.add<std::string>("singleRPReconstructionLabel", "singleRP")
166  ->setComment("output label for single-RP reconstruction products");
167 
168  desc.add<std::string>("multiRPReconstructionLabel", "multiRP")
169  ->setComment("output label for multi-RP reconstruction products");
170 
171  desc.add<double>("localAngleXMin", -0.03)->setComment("minimal accepted value of local horizontal angle (rad)");
172  desc.add<double>("localAngleXMax", +0.03)->setComment("maximal accepted value of local horizontal angle (rad)");
173  desc.add<double>("localAngleYMin", -0.04)->setComment("minimal accepted value of local vertical angle (rad)");
174  desc.add<double>("localAngleYMax", +0.04)->setComment("maximal accepted value of local vertical angle (rad)");
175 
176  std::vector<edm::ParameterSet> config;
177 
178  desc.add<unsigned int>("max_n_timing_tracks", 5)->setComment("maximum number of timing tracks per RP");
179 
180  desc.add<double>("default_time", 0.)->setComment("proton time to be used when no timing information available");
181 
182  desc.add<bool>("fitVtxY", true)
183  ->setComment("for multi-RP reconstruction, flag whether y* should be free fit parameter");
184 
185  desc.add<bool>("useImprovedInitialEstimate", true)
186  ->setComment(
187  "for multi-RP reconstruction, flag whether a quadratic estimate of the initial point should be used");
188 
189  desc.add<std::string>("multiRPAlgorithm", "chi2")
190  ->setComment("algorithm for multi-RP reco, options include chi2, newton, anal-iter");
191 
192  descriptions.add("ctppsProtonsDefault", desc);
193 }
194 
195 //----------------------------------------------------------------------------------------------------
196 
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  // get conditions
211 
214 
216 
218 
219  // re-initialise algorithm upon crossing-angle change
220  if (opticsWatcher_.check(iSetup)) {
221  if (hOpticalFunctions->empty()) {
222  edm::LogInfo("CTPPSProtonProducer") << "No optical functions available, reconstruction disabled.";
224  opticsValid_ = false;
225  } else {
226  algorithm_.init(*hOpticalFunctions);
227  opticsValid_ = true;
228  }
229  }
230 
231  // do reconstruction only if optics is valid
232  if (opticsValid_) {
233  // prepare log
234  std::ostringstream ssLog;
235  if (verbosity_)
236  ssLog << "* input tracks:" << std::endl;
237 
238  // select tracks with small local angles, split them by LHC sector and tracker/timing RPs
239  std::map<unsigned int, std::vector<unsigned int>> trackingSelection, timingSelection;
240 
241  for (unsigned int idx = 0; idx < hTracks->size(); ++idx) {
242  const auto &tr = hTracks->at(idx);
243 
244  if (tr.tx() < localAngleXMin_ || tr.tx() > localAngleXMax_ || tr.ty() < localAngleYMin_ ||
245  tr.ty() > localAngleYMax_)
246  continue;
247 
249  if (tr.pixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::allShiftedPlanes ||
250  tr.pixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::mixedPlanes)
251  continue;
252  }
253 
254  const CTPPSDetId rpId(tr.rpId());
255 
256  if (verbosity_)
257  ssLog << "\t"
258  << "[" << idx << "] " << tr.rpId() << " (" << (rpId.arm() * 100 + rpId.station() * 10 + rpId.rp())
259  << "): "
260  << "x=" << tr.x() << " +- " << tr.xUnc() << " mm, "
261  << "y=" << tr.y() << " +- " << tr.yUnc() << " mm" << std::endl;
262 
263  const bool trackerRP =
264  (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
265 
266  if (trackerRP)
267  trackingSelection[rpId.arm()].push_back(idx);
268  else
269  timingSelection[rpId.arm()].push_back(idx);
270  }
271 
272  // process each arm
273  for (const auto &arm_it : trackingSelection) {
274  const auto &indices = arm_it.second;
275 
276  const auto &ac = ppsAssociationCuts->getAssociationCuts(arm_it.first);
277 
278  // do single-RP reco if needed
279  std::map<unsigned int, reco::ForwardProton> singleRPResultsIndexed;
280  if (doSingleRPReconstruction_ || ac.isApplied(ac.qXi) || ac.isApplied(ac.qThetaY)) {
281  for (const auto &idx : indices) {
282  if (verbosity_)
283  ssLog << std::endl << "* reconstruction from track " << idx << std::endl;
284 
285  singleRPResultsIndexed[idx] =
286  algorithm_.reconstructFromSingleRP(CTPPSLocalTrackLiteRef(hTracks, idx), lhcInfoCombined.energy, 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).rpId());
296 
297  // do multi-RP reco if chosen
298  if (doMultiRPReconstruction_ && rpIds.size() == 2) {
299  // find matching track pairs from different tracking RPs, ordered: i=near, j=far RP
300  std::vector<std::pair<unsigned int, unsigned int>> idx_pairs;
301  std::map<unsigned int, unsigned int> idx_pair_multiplicity;
302  for (const auto &i : indices) {
303  for (const auto &j : indices) {
304  const auto &tr_i = hTracks->at(i);
305  const auto &tr_j = hTracks->at(j);
306 
307  const double z_i = hGeometry->rpTranslation(tr_i.rpId()).z();
308  const double z_j = hGeometry->rpTranslation(tr_j.rpId()).z();
309 
310  const auto &pr_i = singleRPResultsIndexed[i];
311  const auto &pr_j = singleRPResultsIndexed[j];
312 
313  if (tr_i.rpId() == tr_j.rpId())
314  continue;
315 
316  if (std::abs(z_i) >= std::abs(z_j))
317  continue;
318 
319  bool matching = true;
320 
321  if (!ac.isSatisfied(ac.qX, tr_i.x(), tr_i.y(), lhcInfoCombined.crossingAngle(), tr_i.x() - tr_j.x()))
322  matching = false;
323  else if (!ac.isSatisfied(ac.qY, tr_i.x(), tr_i.y(), lhcInfoCombined.crossingAngle(), tr_i.y() - tr_j.y()))
324  matching = false;
325  else if (!ac.isSatisfied(
326  ac.qXi, tr_i.x(), tr_i.y(), lhcInfoCombined.crossingAngle(), pr_i.xi() - pr_j.xi()))
327  matching = false;
328  else if (!ac.isSatisfied(ac.qThetaY,
329  tr_i.x(),
330  tr_i.y(),
331  lhcInfoCombined.crossingAngle(),
332  pr_i.thetaY() - pr_j.thetaY()))
333  matching = false;
334 
335  if (!matching)
336  continue;
337 
338  idx_pairs.emplace_back(i, j);
339  idx_pair_multiplicity[i]++;
340  idx_pair_multiplicity[j]++;
341  }
342  }
343 
344  // evaluate track multiplicity in each timing RP
345  std::map<unsigned int, unsigned int> timing_RP_track_multiplicity;
346  for (const auto &ti : timingSelection[arm_it.first]) {
347  const auto &tr = hTracks->at(ti);
348  timing_RP_track_multiplicity[tr.rpId()]++;
349  }
350 
351  // associate tracking-RP pairs with timing-RP tracks
352  std::map<unsigned int, std::vector<unsigned int>> matched_timing_track_indices;
353  std::map<unsigned int, unsigned int> matched_timing_track_multiplicity;
354  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx) {
355  const auto &i = idx_pairs[pr_idx].first;
356  const auto &j = idx_pairs[pr_idx].second;
357 
358  // skip non-unique associations
359  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
360  continue;
361 
362  const auto &tr_i = hTracks->at(i);
363  const auto &tr_j = hTracks->at(j);
364 
365  const double z_i = hGeometry->rpTranslation(tr_i.rpId()).z();
366  const double z_j = hGeometry->rpTranslation(tr_j.rpId()).z();
367 
368  for (const auto &ti : timingSelection[arm_it.first]) {
369  const auto &tr_ti = hTracks->at(ti);
370 
371  // skip if timing RP saturated (high track multiplicity)
372  if (timing_RP_track_multiplicity[tr_ti.rpId()] > max_n_timing_tracks_)
373  continue;
374 
375  // interpolation from tracking RPs
376  const double z_ti = hGeometry->rpTranslation(tr_ti.rpId()).z();
377  const double f_i = (z_ti - z_j) / (z_i - z_j), f_j = (z_i - z_ti) / (z_i - z_j);
378  const double x_inter = f_i * tr_i.x() + f_j * tr_j.x();
379  const double x_inter_unc_sq =
380  f_i * f_i * tr_i.xUnc() * tr_i.xUnc() + f_j * f_j * tr_j.xUnc() * tr_j.xUnc();
381 
382  const double de_x = tr_ti.x() - x_inter;
383  const double de_x_unc = sqrt(tr_ti.xUnc() * tr_ti.xUnc() + x_inter_unc_sq);
384  const double r = (de_x_unc > 0.) ? de_x / de_x_unc : 1E100;
385 
386  const bool matching = (ac.getTiTrMin() < r && r < ac.getTiTrMax());
387 
388  if (verbosity_)
389  ssLog << "ti=" << ti << ", i=" << i << ", j=" << j << " | z_ti=" << z_ti << ", z_i=" << z_i
390  << ", z_j=" << z_j << " | x_ti=" << tr_ti.x() << ", x_inter=" << x_inter << ", de_x=" << de_x
391  << ", de_x_unc=" << de_x_unc << ", matching=" << matching << std::endl;
392 
393  if (!matching)
394  continue;
395 
396  matched_timing_track_indices[pr_idx].push_back(ti);
397  matched_timing_track_multiplicity[ti]++;
398  }
399  }
400 
401  // process associated tracks
402  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx) {
403  const auto &i = idx_pairs[pr_idx].first;
404  const auto &j = idx_pairs[pr_idx].second;
405 
406  // skip non-unique associations of tracking-RP tracks
407  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
408  continue;
409 
410  if (verbosity_)
411  ssLog << std::endl
412  << "* reconstruction from tracking-RP tracks: " << i << ", " << j << " and timing-RP tracks: ";
413 
414  // buffer contributing tracks
415  CTPPSLocalTrackLiteRefVector sel_tracks;
416  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, i));
417  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, j));
418 
419  CTPPSLocalTrackLiteRefVector sel_track_for_kin_reco = sel_tracks;
420 
421  // process timing-RP data
422  double sw = 0., swt = 0.;
423  for (const auto &ti : matched_timing_track_indices[pr_idx]) {
424  // skip non-unique associations of timing-RP tracks
425  if (matched_timing_track_multiplicity[ti] > 1)
426  continue;
427 
428  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, ti));
429 
430  if (verbosity_)
431  ssLog << ti << ", ";
432 
433  const auto &tr = hTracks->at(ti);
434  const double t_unc = tr.timeUnc();
435  const double w = (t_unc > 0.) ? 1. / t_unc / t_unc : 1.;
436  sw += w;
437  swt += w * tr.time();
438  }
439 
440  float time = default_time_, time_unc = 0.;
441  if (sw > 0.) {
442  time = swt / sw;
443  time_unc = 1. / sqrt(sw);
444  }
445 
446  if (verbosity_)
447  ssLog << std::endl << " time = " << time << " +- " << time_unc << std::endl;
448 
449  // process tracking-RP data
450  reco::ForwardProton proton =
451  algorithm_.reconstructFromMultiRP(sel_track_for_kin_reco, lhcInfoCombined.energy, ssLog);
452 
453  // save combined output
454  proton.setContributingLocalTracks(sel_tracks);
455  proton.setTime(time);
456  proton.setTimeError(time_unc);
457 
458  pOutMultiRP->emplace_back(proton);
459  }
460  }
461 
462  // save single-RP results (un-indexed)
463  for (const auto &p : singleRPResultsIndexed)
464  pOutSingleRP->emplace_back(p.second);
465  }
466 
467  // dump log
468  if (verbosity_)
469  edm::LogInfo("CTPPSProtonProducer") << ssLog.str();
470  }
471  }
472 
473  // save output
475  iEvent.put(std::move(pOutSingleRP), singleRPReconstructionLabel_);
476 
478  iEvent.put(std::move(pOutMultiRP), multiRPReconstructionLabel_);
479 }
480 
481 //----------------------------------------------------------------------------------------------------
482 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::ESGetToken< LHCInfoPerFill, LHCInfoPerFillRcd > lhcInfoPerFillToken_
reco::ForwardProton reconstructFromSingleRP(const CTPPSLocalTrackLiteRef &track, const float energy, std::ostream &os) const
run proton reconstruction using single-RP strategy
T w() const
const std::string opticsLabel_
float crossingAngle() const
const std::string ppsAssociationCutsLabel_
~CTPPSProtonProducer() override=default
const edm::ESGetToken< PPSAssociationCuts, PPSAssociationCutsRcd > ppsAssociationCutsToken_
unsigned int max_n_timing_tracks_
Event setup record containing the real (actual) geometry information.
const edm::ESGetToken< LHCInfoPerLS, LHCInfoPerLSRcd > lhcInfoPerLSToken_
void init(const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions)
const CutsPerArm & getAssociationCuts(const int sector) const
int iEvent
Definition: GenABIO.cc:224
const std::string lhcInfoPerLSLabel_
void setTimeError(float time_err)
const std::string lhcInfoLabel_
const edm::ESGetToken< LHCInfo, LHCInfoRcd > lhcInfoToken_
T sqrt(T t)
Definition: SSEVec.h:23
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.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setTime(float time)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
reco::ForwardProton reconstructFromMultiRP(const CTPPSLocalTrackLiteRefVector &tracks, const float energy, std::ostream &os) const
run proton reconstruction using multiple-RP strategy
The manager class for TOTEM RP geometry.
Definition: CTPPSGeometry.h:30
Log< level::Info, false > LogInfo
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
Vector rpTranslation(unsigned int id) 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:57
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLT enums.
ProtonReconstructionAlgorithm algorithm_
dictionary config
Read in AllInOne config in JSON format.
Definition: DiMuonV_cfg.py:30
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
const std::string lhcInfoPerFillLabel_
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
const edm::ESGetToken< CTPPSGeometry, VeryForwardRealGeometryRecord > geometryToken_
std::string singleRPReconstructionLabel_
def move(src, dest)
Definition: eostools.py:511
const edm::ESGetToken< LHCInterpolatedOpticalFunctionsSetCollection, CTPPSInterpolatedOpticsRcd > opticalFunctionsToken_
std::string multiRPReconstructionLabel_
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tracksToken_