CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
CTPPSProtonProducer Class Reference
Inheritance diagram for CTPPSProtonProducer:
edm::stream::EDProducer<>

Classes

struct  AssociationCuts
 

Public Member Functions

 CTPPSProtonProducer (const edm::ParameterSet &)
 
 ~CTPPSProtonProducer () override=default
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

ProtonReconstructionAlgorithm algorithm_
 
std::map< unsigned int, AssociationCutsassociation_cuts_
 
bool doMultiRPReconstruction_
 
bool doSingleRPReconstruction_
 
std::string lhcInfoLabel_
 
double localAngleXMax_
 
double localAngleXMin_
 
double localAngleYMax_
 
double localAngleYMin_
 
unsigned int max_n_timing_tracks_
 
std::string multiRPReconstructionLabel_
 
std::string opticsLabel_
 
bool opticsValid_
 
edm::ESWatcher< CTPPSInterpolatedOpticsRcdopticsWatcher_
 
std::string singleRPReconstructionLabel_
 
edm::EDGetTokenT< CTPPSLocalTrackLiteCollectiontracksToken_
 
unsigned int verbosity_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Definition at line 35 of file CTPPSProtonProducer.cc.

Constructor & Destructor Documentation

CTPPSProtonProducer::CTPPSProtonProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 122 of file CTPPSProtonProducer.cc.

References association_cuts_, doMultiRPReconstruction_, doSingleRPReconstruction_, edm::ParameterSet::getParameterSet(), multiRPReconstructionLabel_, singleRPReconstructionLabel_, and AlCaHLTBitMon_QueryRunRegistry::string.

122  :
123  tracksToken_ (consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagLocalTrackLite"))),
124  lhcInfoLabel_ (iConfig.getParameter<std::string>("lhcInfoLabel")),
125  opticsLabel_ (iConfig.getParameter<std::string>("opticsLabel")),
126  verbosity_ (iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)),
127  doSingleRPReconstruction_ (iConfig.getParameter<bool>("doSingleRPReconstruction")),
128  doMultiRPReconstruction_ (iConfig.getParameter<bool>("doMultiRPReconstruction")),
129  singleRPReconstructionLabel_(iConfig.getParameter<std::string>("singleRPReconstructionLabel")),
130  multiRPReconstructionLabel_ (iConfig.getParameter<std::string>("multiRPReconstructionLabel")),
131 
132  localAngleXMin_ (iConfig.getParameter<double>("localAngleXMin")),
133  localAngleXMax_ (iConfig.getParameter<double>("localAngleXMax")),
134  localAngleYMin_ (iConfig.getParameter<double>("localAngleYMin")),
135  localAngleYMax_ (iConfig.getParameter<double>("localAngleYMax")),
136 
137  max_n_timing_tracks_ (iConfig.getParameter<unsigned int>("max_n_timing_tracks")),
138 
139  algorithm_ (iConfig.getParameter<bool>("fitVtxY"), iConfig.getParameter<bool>("useImprovedInitialEstimate"), verbosity_),
140  opticsValid_(false)
141 {
142  for (const std::string &sector : { "45", "56" })
143  {
144  const unsigned int arm = (sector == "45") ? 0 : 1;
145  association_cuts_[arm].load(iConfig.getParameterSet("association_cuts_" + sector));
146  }
147 
149  produces<reco::ForwardProtonCollection>(singleRPReconstructionLabel_);
150 
152  produces<reco::ForwardProtonCollection>(multiRPReconstructionLabel_);
153 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
unsigned int max_n_timing_tracks_
std::map< unsigned int, AssociationCuts > association_cuts_
ParameterSet const & getParameterSet(std::string const &) const
ProtonReconstructionAlgorithm algorithm_
std::string singleRPReconstructionLabel_
std::string multiRPReconstructionLabel_
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tracksToken_
CTPPSProtonProducer::~CTPPSProtonProducer ( )
overridedefault

Member Function Documentation

void CTPPSProtonProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 157 of file CTPPSProtonProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), looper::config, CTPPSProtonProducer::AssociationCuts::getDefaultParameters(), edm::ParameterDescriptionNode::setComment(), and AlCaHLTBitMon_QueryRunRegistry::string.

158 {
160 
161  desc.add<edm::InputTag>("tagLocalTrackLite", edm::InputTag("ctppsLocalTrackLiteProducer"))
162  ->setComment("specification of the input lite-track collection");
163 
164  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
165  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics record");
166 
167  desc.addUntracked<unsigned int>("verbosity", 0)->setComment("verbosity level");
168 
169  desc.add<bool>("doSingleRPReconstruction", true)
170  ->setComment("flag whether to apply single-RP reconstruction strategy");
171 
172  desc.add<bool>("doMultiRPReconstruction", true)
173  ->setComment("flag whether to apply multi-RP reconstruction strategy");
174 
175  desc.add<std::string>("singleRPReconstructionLabel", "singleRP")
176  ->setComment("output label for single-RP reconstruction products");
177 
178  desc.add<std::string>("multiRPReconstructionLabel", "multiRP")
179  ->setComment("output label for multi-RP reconstruction products");
180 
181  desc.add<double>("localAngleXMin", -0.03)->setComment("minimal accepted value of local horizontal angle (rad)");
182  desc.add<double>("localAngleXMax", +0.03)->setComment("maximal accepted value of local horizontal angle (rad)");
183  desc.add<double>("localAngleYMin", -0.04)->setComment("minimal accepted value of local vertical angle (rad)");
184  desc.add<double>("localAngleYMax", +0.04)->setComment("maximal accepted value of local vertical angle (rad)");
185 
186  for (const std::string &sector : { "45", "56" })
187  {
188  desc.add<edm::ParameterSetDescription>("association_cuts_" + sector, AssociationCuts::getDefaultParameters())
189  ->setComment("track-association cuts for sector " + sector);
190  }
191 
192  std::vector<edm::ParameterSet> config;
193 
194  desc.add<unsigned int>("max_n_timing_tracks", 5)->setComment("maximum number of timing tracks per RP");
195 
196  desc.add<bool>("fitVtxY", true)
197  ->setComment("for multi-RP reconstruction, flag whether y* should be free fit parameter");
198 
199  desc.add<bool>("useImprovedInitialEstimate", true)
200  ->setComment("for multi-RP reconstruction, flag whether a quadratic estimate of the initial point should be used");
201 
202  descriptions.add("ctppsProtons", desc);
203 }
void setComment(std::string const &value)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
config
Definition: looper.py:291
static edm::ParameterSetDescription getDefaultParameters()
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void CTPPSProtonProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 207 of file CTPPSProtonProducer.cc.

References funct::abs(), algorithm_, association_cuts_, edm::ESWatcher< T >::check(), alignCSCRings::de_x, DEFINE_FWK_MODULE, doMultiRPReconstruction_, doSingleRPReconstruction_, edm::EventSetup::get(), edm::Event::getByToken(), CTPPSGeometry::getRPTranslation(), mps_fire::i, training_settings::idx, ProtonReconstructionAlgorithm::init(), lhcInfoLabel_, localAngleXMax_, localAngleXMin_, localAngleYMax_, localAngleYMin_, max_n_timing_tracks_, eostools::move(), multiRPReconstructionLabel_, opticsLabel_, opticsValid_, opticsWatcher_, AlCaHLTBitMon_ParallelJobs::p, edm::RefVector< C, T, F >::push_back(), edm::Event::put(), alignCSCRings::r, ProtonReconstructionAlgorithm::reconstructFromMultiRP(), ProtonReconstructionAlgorithm::reconstructFromSingleRP(), ProtonReconstructionAlgorithm::release(), year_2016_postTS2_cff::rpId, year_2016_cff::rpIds, CTPPSDetId::sdTrackingPixel, CTPPSDetId::sdTrackingStrip, reco::ForwardProton::setContributingLocalTracks(), reco::ForwardProton::setTime(), reco::ForwardProton::setTimeError(), singleRPReconstructionLabel_, mathSSE::sqrt(), ntuplemaker::time, tracksToken_, verbosity_, w, and z.

208 {
209  // get input
211  iEvent.getByToken(tracksToken_, hTracks);
212 
213  // book output
214  std::unique_ptr<reco::ForwardProtonCollection> pOutSingleRP(new reco::ForwardProtonCollection);
215  std::unique_ptr<reco::ForwardProtonCollection> pOutMultiRP(new reco::ForwardProtonCollection);
216 
217  // continue only if there is something to process
218  // NB: this avoids loading (possibly non-existing) conditions in workflows without proton data
219  if (!hTracks->empty())
220  {
221  // get conditions
222  edm::ESHandle<LHCInfo> hLHCInfo;
223  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
224 
226  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(opticsLabel_, hOpticalFunctions);
227 
229  iSetup.get<VeryForwardRealGeometryRecord>().get(hGeometry);
230 
231  // re-initialise algorithm upon crossing-angle change
232  if (opticsWatcher_.check(iSetup))
233  {
234  if (hOpticalFunctions->empty()) {
235  edm::LogInfo("CTPPSProtonProducer") << "No optical functions available, reconstruction disabled.";
237  opticsValid_ = false;
238  }
239  else {
240  algorithm_.init(*hOpticalFunctions);
241  opticsValid_ = true;
242  }
243  }
244 
245  // do reconstruction only if optics is valid
246  if (opticsValid_)
247  {
248  // prepare log
249  std::ostringstream ssLog;
250  if (verbosity_)
251  ssLog << "* input tracks:" << std::endl;
252 
253  // select tracks with small local angles, split them by LHC sector and tracker/timing RPs
254  std::map<unsigned int, std::vector<unsigned int>> trackingSelection, timingSelection;
255 
256  for (unsigned int idx = 0; idx < hTracks->size(); ++idx)
257  {
258  const auto& tr = hTracks->at(idx);
259 
260  if (tr.getTx() < localAngleXMin_ || tr.getTx() > localAngleXMax_
261  || tr.getTy() < localAngleYMin_ || tr.getTy() > localAngleYMax_)
262  continue;
263 
264  const CTPPSDetId rpId(tr.getRPId());
265 
266  if (verbosity_)
267  ssLog << "\t"
268  << "[" << idx << "] "
269  << tr.getRPId() << " (" << (rpId.arm()*100 + rpId.station()*10 + rpId.rp()) << "): "
270  << "x=" << tr.getX() << " +- " << tr.getXUnc() << " mm, "
271  << "y=" << tr.getY() << " +- " << tr.getYUnc() << " mm" << std::endl;
272 
273  const bool trackerRP = (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
274 
275  if (trackerRP)
276  trackingSelection[rpId.arm()].push_back(idx);
277  else
278  timingSelection[rpId.arm()].push_back(idx);
279  }
280 
281  // process each arm
282  for (const auto &arm_it : trackingSelection)
283  {
284  const auto &indices = arm_it.second;
285 
286  const auto &ac = association_cuts_[arm_it.first];
287 
288  // do single-RP reco if needed
289  std::map<unsigned int, reco::ForwardProton> singleRPResultsIndexed;
290  if (doSingleRPReconstruction_ || ac.xi_cut_apply || ac.th_y_cut_apply)
291  {
292  for (const auto &idx : indices)
293  {
294  if (verbosity_)
295  ssLog << std::endl << "* reconstruction from track " << idx << std::endl;
296 
297  singleRPResultsIndexed[idx] = algorithm_.reconstructFromSingleRP(CTPPSLocalTrackLiteRef(hTracks, idx), *hLHCInfo, ssLog);
298  }
299  }
300 
301  // check that exactly two tracking RPs are involved
302  // - 1 is insufficient for multi-RP reconstruction
303  // - PPS did not use more than 2 tracking RPs per arm -> algorithms are tuned to this
304  std::set<unsigned int> rpIds;
305  for (const auto &idx : indices)
306  rpIds.insert(hTracks->at(idx).getRPId());
307 
308  // do multi-RP reco if chosen
309  if (doMultiRPReconstruction_ && rpIds.size() == 2)
310  {
311  // find matching track pairs from different tracking RPs
312  std::vector<std::pair<unsigned int, unsigned int>> idx_pairs;
313  std::map<unsigned int, unsigned int> idx_pair_multiplicity;
314  for (const auto &i : indices)
315  {
316  for (const auto &j : indices)
317  {
318  if (j <= i)
319  continue;
320 
321  const auto &tr_i = hTracks->at(i);
322  const auto &tr_j = hTracks->at(j);
323 
324  const auto &pr_i = singleRPResultsIndexed[i];
325  const auto &pr_j = singleRPResultsIndexed[j];
326 
327  if (tr_i.getRPId() == tr_j.getRPId())
328  continue;
329 
330  bool matching = true;
331 
332  if (ac.x_cut_apply && std::abs(tr_i.getX() - tr_j.getX()) > ac.x_cut_value)
333  matching = false;
334  if (ac.y_cut_apply && std::abs(tr_i.getY() - tr_j.getY()) > ac.y_cut_value)
335  matching = false;
336  if (ac.xi_cut_apply && std::abs(pr_i.xi() - pr_j.xi()) > ac.xi_cut_value)
337  matching = false;
338  if (ac.th_y_cut_apply && std::abs(pr_i.thetaY() - pr_j.thetaY()) > ac.th_y_cut_value)
339  matching = false;
340 
341  if (!matching)
342  continue;
343 
344  idx_pairs.emplace_back(i, j);
345  idx_pair_multiplicity[i]++;
346  idx_pair_multiplicity[j]++;
347  }
348  }
349 
350  // evaluate track multiplicity in each timing RP
351  std::map<unsigned int, unsigned int> timing_RP_track_multiplicity;
352  for (const auto &ti : timingSelection[arm_it.first])
353  {
354  const auto &tr = hTracks->at(ti);
355  timing_RP_track_multiplicity[tr.getRPId()]++;
356  }
357 
358  // associate tracking-RP pairs with timing-RP tracks
359  std::map<unsigned int, std::vector<unsigned int>> matched_timing_track_indices;
360  std::map<unsigned int, unsigned int> matched_timing_track_multiplicity;
361  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx)
362  {
363  const auto &i = idx_pairs[pr_idx].first;
364  const auto &j = idx_pairs[pr_idx].second;
365 
366  // skip non-unique associations
367  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
368  continue;
369 
370  const auto &tr_i = hTracks->at(i);
371  const auto &tr_j = hTracks->at(j);
372 
373  const double z_i = hGeometry->getRPTranslation(tr_i.getRPId()).z();
374  const double z_j = hGeometry->getRPTranslation(tr_j.getRPId()).z();
375 
376  for (const auto &ti : timingSelection[arm_it.first])
377  {
378  const auto &tr_ti = hTracks->at(ti);
379 
380  // skip if timing RP saturated (high track multiplicity)
381  if (timing_RP_track_multiplicity[tr_ti.getRPId()] > max_n_timing_tracks_)
382  continue;
383 
384  // interpolation from tracking RPs
385  const double z_ti = - hGeometry->getRPTranslation(tr_ti.getRPId()).z(); // the minus sign fixes a bug in the diamond geometry
386  const double f_i = (z_ti - z_j) / (z_i - z_j), f_j = (z_i - z_ti) / (z_i - z_j);
387  const double x_inter = f_i * tr_i.getX() + f_j * tr_j.getX();
388  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();
389 
390  const double de_x = tr_ti.getX() - x_inter;
391  const double de_x_unc = sqrt(tr_ti.getXUnc()*tr_ti.getXUnc() + x_inter_unc_sq);
392  const double r = (de_x_unc > 0.) ? de_x / de_x_unc : 1E100;
393 
394  const bool matching = (ac.ti_tr_min < r && r < ac.ti_tr_max);
395 
396  if (verbosity_)
397  ssLog << "ti=" << ti << ", i=" << i << ", j=" << j
398  << " | z_ti=" << z_ti << ", z_i=" << z_i << ", z_j=" << z_j
399  << " | x_ti=" << tr_ti.getX() << ", x_inter=" << x_inter << ", de_x=" << de_x << ", de_x_unc=" << de_x_unc
400  << ", matching=" << matching << std::endl;
401 
402  if (!matching)
403  continue;
404 
405  matched_timing_track_indices[pr_idx].push_back(ti);
406  matched_timing_track_multiplicity[ti]++;
407  }
408  }
409 
410  // process associated tracks
411  for (unsigned int pr_idx = 0; pr_idx < idx_pairs.size(); ++pr_idx)
412  {
413  const auto &i = idx_pairs[pr_idx].first;
414  const auto &j = idx_pairs[pr_idx].second;
415 
416  // skip non-unique associations of tracking-RP tracks
417  if (idx_pair_multiplicity[i] > 1 || idx_pair_multiplicity[j] > 1)
418  continue;
419 
420  if (verbosity_)
421  ssLog << std::endl << "* reconstruction from tracking-RP tracks: " << i << ", " << j << " and timing-RP tracks: ";
422 
423  // buffer contributing tracks
424  CTPPSLocalTrackLiteRefVector sel_tracks;
425  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, i));
426  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, j));
427 
428  CTPPSLocalTrackLiteRefVector sel_track_for_kin_reco = sel_tracks;
429 
430  // process timing-RP data
431  double sw=0., swt=0.;
432  for (const auto &ti : matched_timing_track_indices[pr_idx])
433  {
434  // skip non-unique associations of timing-RP tracks
435  if (matched_timing_track_multiplicity[ti] > 1)
436  continue;
437 
438  sel_tracks.push_back(CTPPSLocalTrackLiteRef(hTracks, ti));
439 
440  if (verbosity_)
441  ssLog << ti << ", ";
442 
443  const auto &tr = hTracks->at(ti);
444  const double t_unc = tr.getTimeUnc();
445  const double w = (t_unc > 0.) ? 1./t_unc/t_unc : 1.;
446  sw += w;
447  swt += w * tr.getTime();
448  }
449 
450  float time = 0., time_unc = 0.;
451  if (sw > 0.)
452  {
453  time = swt / sw;
454  time_unc = 1. / sqrt(sw);
455  }
456 
457  if (verbosity_)
458  ssLog << std::endl << " time = " << time << " +- " << time_unc << std::endl;
459 
460  // process tracking-RP data
461  reco::ForwardProton proton = algorithm_.reconstructFromMultiRP(sel_track_for_kin_reco, *hLHCInfo, ssLog);
462 
463  // save combined output
464  proton.setContributingLocalTracks(sel_tracks);
465  proton.setTime(time);
466  proton.setTimeError(time_unc);
467 
468  pOutMultiRP->emplace_back(proton);
469  }
470  }
471 
472  // save single-RP results (un-indexed)
473  for (const auto &p : singleRPResultsIndexed)
474  pOutSingleRP->emplace_back(std::move(p.second));
475  }
476 
477  // dump log
478  if (verbosity_)
479  edm::LogInfo("CTPPSProtonProducer") << ssLog.str();
480  }
481  }
482 
483  // save output
485  iEvent.put(std::move(pOutSingleRP), singleRPReconstructionLabel_);
486 
488  iEvent.put(std::move(pOutMultiRP), multiRPReconstructionLabel_);
489 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
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
unsigned int max_n_timing_tracks_
Event setup record containing the real (actual) geometry information.
void init(const LHCInterpolatedOpticalFunctionsSetCollection &opticalFunctions)
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
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)
CLHEP::Hep3Vector getRPTranslation(unsigned int id) const
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
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_

Member Data Documentation

ProtonReconstructionAlgorithm CTPPSProtonProducer::algorithm_
private

Definition at line 114 of file CTPPSProtonProducer.cc.

Referenced by produce().

std::map<unsigned int, AssociationCuts> CTPPSProtonProducer::association_cuts_
private

Definition at line 110 of file CTPPSProtonProducer.cc.

Referenced by CTPPSProtonProducer(), and produce().

bool CTPPSProtonProducer::doMultiRPReconstruction_
private

Definition at line 54 of file CTPPSProtonProducer.cc.

Referenced by CTPPSProtonProducer(), and produce().

bool CTPPSProtonProducer::doSingleRPReconstruction_
private

Definition at line 53 of file CTPPSProtonProducer.cc.

Referenced by CTPPSProtonProducer(), and produce().

std::string CTPPSProtonProducer::lhcInfoLabel_
private

Definition at line 48 of file CTPPSProtonProducer.cc.

Referenced by produce().

double CTPPSProtonProducer::localAngleXMax_
private

Definition at line 59 of file CTPPSProtonProducer.cc.

Referenced by produce().

double CTPPSProtonProducer::localAngleXMin_
private

Definition at line 59 of file CTPPSProtonProducer.cc.

Referenced by produce().

double CTPPSProtonProducer::localAngleYMax_
private

Definition at line 59 of file CTPPSProtonProducer.cc.

Referenced by produce().

double CTPPSProtonProducer::localAngleYMin_
private

Definition at line 59 of file CTPPSProtonProducer.cc.

Referenced by produce().

unsigned int CTPPSProtonProducer::max_n_timing_tracks_
private

Definition at line 112 of file CTPPSProtonProducer.cc.

Referenced by produce().

std::string CTPPSProtonProducer::multiRPReconstructionLabel_
private

Definition at line 57 of file CTPPSProtonProducer.cc.

Referenced by CTPPSProtonProducer(), and produce().

std::string CTPPSProtonProducer::opticsLabel_
private

Definition at line 49 of file CTPPSProtonProducer.cc.

Referenced by produce().

bool CTPPSProtonProducer::opticsValid_
private

Definition at line 116 of file CTPPSProtonProducer.cc.

Referenced by produce().

edm::ESWatcher<CTPPSInterpolatedOpticsRcd> CTPPSProtonProducer::opticsWatcher_
private

Definition at line 117 of file CTPPSProtonProducer.cc.

Referenced by produce().

std::string CTPPSProtonProducer::singleRPReconstructionLabel_
private

Definition at line 56 of file CTPPSProtonProducer.cc.

Referenced by CTPPSProtonProducer(), and produce().

edm::EDGetTokenT<CTPPSLocalTrackLiteCollection> CTPPSProtonProducer::tracksToken_
private

Definition at line 46 of file CTPPSProtonProducer.cc.

Referenced by produce().

unsigned int CTPPSProtonProducer::verbosity_
private

Definition at line 51 of file CTPPSProtonProducer.cc.

Referenced by produce().