CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Public Attributes
PPSAlignmentWorker::SectorData Struct Reference

Classes

struct  SlicePlots
 

Public Member Functions

void init (DQMStore::IBooker &iBooker, const PPSAlignmentConfiguration &cfg, const PPSAlignmentConfiguration::SectorConfig &scfg, const std::string &rootDir, bool debug)
 
unsigned int process (const CTPPSLocalTrackLiteCollection &tracks, const PPSAlignmentConfiguration &cfg, bool debug)
 

Public Attributes

MonitorElementh2_cut_h_aft
 
MonitorElementh2_cut_h_bef
 
MonitorElementh2_cut_v_aft
 
MonitorElementh2_cut_v_bef
 
MonitorElementh_q_cut_h_aft
 
MonitorElementh_q_cut_h_bef
 
MonitorElementh_q_cut_v_aft
 
MonitorElementh_q_cut_v_bef
 
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_aft_sel
 
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_bef_sel
 
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_mlt_sel
 
MonitorElementp_x_diffFN_vs_x_N
 
MonitorElementp_y_diffFN_vs_y_F
 
PPSAlignmentConfiguration::SectorConfig scfg_
 
std::map< unsigned int, SlicePlotsx_slice_plots_F
 
std::map< unsigned int, SlicePlotsx_slice_plots_N
 

Detailed Description

Definition at line 46 of file PPSAlignmentWorker.cc.

Member Function Documentation

◆ init()

void PPSAlignmentWorker::SectorData::init ( DQMStore::IBooker iBooker,
const PPSAlignmentConfiguration cfg,
const PPSAlignmentConfiguration::SectorConfig scfg,
const std::string &  rootDir,
bool  debug 
)

Definition at line 203 of file PPSAlignmentWorker.cc.

References dqm::implementation::IBooker::book1DD(), dqm::implementation::IBooker::book2DD(), dqm::implementation::IBooker::bookProfile(), visDQMUpload::buf, looper::cfg, debug, h2_cut_h_aft, h2_cut_h_bef, h2_cut_v_aft, h2_cut_v_bef, h_q_cut_h_aft, h_q_cut_h_bef, h_q_cut_v_aft, h_q_cut_v_bef, mps_fire::i, PPSAlignmentConfiguration::RPConfig::id_, m_h2_y_vs_x_aft_sel, m_h2_y_vs_x_bef_sel, m_h2_y_vs_x_mlt_sel, PPSAlignmentConfiguration::RPConfig::name_, PPSAlignmentConfiguration::SectorConfig::name_, p_x_diffFN_vs_x_N, python.rootplot.root2matplotlib::replace(), PPSAlignmentConfiguration::SectorConfig::rp_F_, PPSAlignmentConfiguration::SectorConfig::rp_N_, scfg_, dqm::implementation::NavigatorBase::setCurrentFolder(), PPSAlignmentConfiguration::RPConfig::x_slice_min_, PPSAlignmentConfiguration::RPConfig::x_slice_n_, x_slice_plots_F, x_slice_plots_N, and PPSAlignmentConfiguration::RPConfig::x_slice_w_.

Referenced by PPSAlignmentWorker::bookHistograms().

207  {
208  scfg_ = scfg;
209 
210  // binning
211  const double bin_size_x = cfg.binning().bin_size_x_;
212  const unsigned int n_bins_x = cfg.binning().n_bins_x_;
213 
214  const double pixel_x_offset = cfg.binning().pixel_x_offset_;
215 
216  const double x_min_pix = pixel_x_offset, x_max_pix = pixel_x_offset + n_bins_x * bin_size_x;
217  const double x_min_str = 0., x_max_str = n_bins_x * bin_size_x;
218 
219  const unsigned int n_bins_y = cfg.binning().n_bins_y_;
220  const double y_min = cfg.binning().y_min_, y_max = cfg.binning().y_max_;
221 
222  // hit distributions
223  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/before selection/" + scfg_.rp_N_.name_);
225  iBooker.book2DD("h2_y_vs_x", ";x;y", n_bins_x, x_min_str, x_max_str, n_bins_y, y_min, y_max);
226  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/before selection/" + scfg_.rp_F_.name_);
228  iBooker.book2DD("h2_y_vs_x", ";x;y", n_bins_x, x_min_pix, x_max_pix, n_bins_y, y_min, y_max);
229 
230  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/multiplicity selection/" + scfg_.rp_N_.name_);
232  iBooker.book2DD("h2_y_vs_x", ";x;y", n_bins_x, x_min_str, x_max_str, n_bins_y, y_min, y_max);
233  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/multiplicity selection/" + scfg_.rp_F_.name_);
235  iBooker.book2DD("h2_y_vs_x", ";x;y", n_bins_x, x_min_pix, x_max_pix, n_bins_y, y_min, y_max);
236 
237  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/after selection/" + scfg_.rp_N_.name_);
239  iBooker.book2DD("h2_y_vs_x", ";x;y", n_bins_x, x_min_str, x_max_str, n_bins_y, y_min, y_max);
240  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/after selection/" + scfg_.rp_F_.name_);
242  iBooker.book2DD("h2_y_vs_x", ";x;y", n_bins_x, x_min_pix, x_max_pix, n_bins_y, y_min, y_max);
243 
244  // cut plots
245  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/cuts/cut_h");
246  h_q_cut_h_bef = iBooker.book1DD("h_q_cut_h_bef", ";cq_h", 400, -2., 2.);
247  h_q_cut_h_aft = iBooker.book1DD("h_q_cut_h_aft", ";cq_h", 400, -2., 2.);
248  h2_cut_h_bef =
249  iBooker.book2DD("h2_cut_h_bef", ";x_up;x_dw", n_bins_x, x_min_str, x_max_str, n_bins_x, x_min_pix, x_max_pix);
250  h2_cut_h_aft =
251  iBooker.book2DD("h2_cut_h_aft", ";x_up;x_dw", n_bins_x, x_min_str, x_max_str, n_bins_x, x_min_pix, x_max_pix);
252 
253  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/cuts/cut_v");
254  h_q_cut_v_bef = iBooker.book1DD("h_q_cut_v_bef", ";cq_v", 400, -2., 2.);
255  h_q_cut_v_aft = iBooker.book1DD("h_q_cut_v_aft", ";cq_v", 400, -2., 2.);
256  h2_cut_v_bef = iBooker.book2DD("h2_cut_v_bef", ";y_up;y_dw", n_bins_y, y_min, y_max, n_bins_y, y_min, y_max);
257  h2_cut_v_aft = iBooker.book2DD("h2_cut_v_aft", ";y_up;y_dw", n_bins_y, y_min, y_max, n_bins_y, y_min, y_max);
258 
259  // near-far plots
260  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/near_far");
261 
262  auto profilePtr = std::make_unique<TProfile>("",
263  ";x_{N};x_{F} - x_{N}",
264  cfg.binning().diffFN_n_bins_x_,
265  cfg.binning().diffFN_x_min_,
266  cfg.binning().diffFN_x_max_);
267  p_x_diffFN_vs_x_N = iBooker.bookProfile("p_x_diffFN_vs_x_N", profilePtr.get());
268 
269  // slice plots
270  for (int i = 0; i < scfg_.rp_N_.x_slice_n_; i++) {
271  const double x_min = scfg_.rp_N_.x_slice_min_ + i * scfg_.rp_N_.x_slice_w_;
272  const double x_max = scfg_.rp_N_.x_slice_min_ + (i + 1) * scfg_.rp_N_.x_slice_w_;
273 
274  char buf[100];
275  sprintf(buf, "%.1f-%.1f", x_min, x_max);
276  std::replace(buf, buf + strlen(buf), '.', '_'); // replace . with _
277  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/near_far/x slices N/" + buf);
278  x_slice_plots_N.insert({i, SlicePlots(iBooker, cfg, debug)});
279  }
280 
281  for (int i = 0; i < scfg_.rp_F_.x_slice_n_; i++) {
282  const double x_min = scfg_.rp_F_.x_slice_min_ + i * scfg_.rp_F_.x_slice_w_;
283  const double x_max = scfg_.rp_F_.x_slice_min_ + (i + 1) * scfg_.rp_F_.x_slice_w_;
284 
285  char buf[100];
286  sprintf(buf, "%.1f-%.1f", x_min, x_max);
287  std::replace(buf, buf + strlen(buf), '.', '_'); // replace . with _
288  iBooker.setCurrentFolder(rootDir + "/" + scfg_.name_ + "/near_far/x slices F/" + buf);
289  x_slice_plots_F.insert({i, SlicePlots(iBooker, cfg, debug)});
290  }
291 }
std::map< unsigned int, SlicePlots > x_slice_plots_F
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_bef_sel
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
def replace(string, replacements)
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_aft_sel
MonitorElement * book1DD(TString const &name, TString const &title, int nchX, double lowX, double highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:155
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
PPSAlignmentConfiguration::SectorConfig scfg_
#define debug
Definition: HDRShower.cc:19
MonitorElement * book2DD(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:338
std::map< unsigned int, SlicePlots > x_slice_plots_N
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_mlt_sel

◆ process()

unsigned int PPSAlignmentWorker::SectorData::process ( const CTPPSLocalTrackLiteCollection tracks,
const PPSAlignmentConfiguration cfg,
bool  debug 
)

Definition at line 293 of file PPSAlignmentWorker.cc.

References CTPPSDetId::arm(), looper::cfg, PPSAlignmentConfiguration::SectorConfig::cut_h_a_, PPSAlignmentConfiguration::SectorConfig::cut_h_apply_, PPSAlignmentConfiguration::SectorConfig::cut_h_c_, PPSAlignmentConfiguration::SectorConfig::cut_h_si_, PPSAlignmentConfiguration::SectorConfig::cut_v_a_, PPSAlignmentConfiguration::SectorConfig::cut_v_apply_, PPSAlignmentConfiguration::SectorConfig::cut_v_c_, PPSAlignmentConfiguration::SectorConfig::cut_v_si_, debug, dqm::impl::MonitorElement::Fill(), h2_cut_h_aft, h2_cut_h_bef, h2_cut_v_aft, h2_cut_v_bef, h_q_cut_h_aft, h_q_cut_h_bef, h_q_cut_v_aft, h_q_cut_v_bef, PPSAlignmentConfiguration::RPConfig::id_, heavyIonCSV_trainingSettings::idx, if(), m_h2_y_vs_x_aft_sel, m_h2_y_vs_x_bef_sel, m_h2_y_vs_x_mlt_sel, p_x_diffFN_vs_x_N, PPSAlignmentConfiguration::SectorConfig::rp_F_, PPSAlignmentConfiguration::SectorConfig::rp_N_, scfg_, tracks, PPSAlignmentConfiguration::RPConfig::x_slice_min_, PPSAlignmentConfiguration::RPConfig::x_slice_n_, x_slice_plots_F, x_slice_plots_N, and PPSAlignmentConfiguration::RPConfig::x_slice_w_.

Referenced by PPSAlignmentWorker::analyze().

295  {
296  CTPPSLocalTrackLiteCollection tracksUp, tracksDw;
297 
298  for (const auto& tr : tracks) {
299  CTPPSDetId rpId(tr.rpId());
300  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
301 
302  if (rpDecId != scfg_.rp_N_.id_ && rpDecId != scfg_.rp_F_.id_)
303  continue;
304 
305  // store the track in the right collection
306  if (rpDecId == scfg_.rp_N_.id_)
307  tracksUp.push_back(tr);
308  if (rpDecId == scfg_.rp_F_.id_)
309  tracksDw.push_back(tr);
310  }
311 
312  // update plots before selection
313  for (const auto& tr : tracksUp)
314  m_h2_y_vs_x_bef_sel[scfg_.rp_N_.id_]->Fill(tr.x(), tr.y());
315 
316  for (const auto& tr : tracksDw)
317  m_h2_y_vs_x_bef_sel[scfg_.rp_F_.id_]->Fill(tr.x(), tr.y());
318 
319  // skip crowded events (multiplicity selection)
320  if (tracksUp.size() < cfg.minRPTracksSize() || tracksUp.size() > cfg.maxRPTracksSize())
321  return 0;
322 
323  if (tracksDw.size() < cfg.minRPTracksSize() || tracksDw.size() > cfg.maxRPTracksSize())
324  return 0;
325 
326  // update plots with multiplicity selection
327  for (const auto& tr : tracksUp)
328  m_h2_y_vs_x_mlt_sel[scfg_.rp_N_.id_]->Fill(tr.x(), tr.y());
329 
330  for (const auto& tr : tracksDw)
331  m_h2_y_vs_x_mlt_sel[scfg_.rp_F_.id_]->Fill(tr.x(), tr.y());
332 
333  // do the selection
334  unsigned int pairsSelected = 0;
335  for (const auto& trUp : tracksUp) {
336  for (const auto& trDw : tracksDw) {
337  h2_cut_h_bef->Fill(trUp.x(), trDw.x());
338  h2_cut_v_bef->Fill(trUp.y(), trDw.y());
339 
340  // horizontal cut
341  const double cq_h = trDw.x() + scfg_.cut_h_a_ * trUp.x() + scfg_.cut_h_c_;
342  h_q_cut_h_bef->Fill(cq_h);
343  const bool cv_h = (std::fabs(cq_h) < cfg.n_si() * scfg_.cut_h_si_);
344 
345  // vertical cut
346  const double cq_v = trDw.y() + scfg_.cut_v_a_ * trUp.y() + scfg_.cut_v_c_;
347  h_q_cut_v_bef->Fill(cq_v);
348  const bool cv_v = (std::fabs(cq_v) < cfg.n_si() * scfg_.cut_v_si_);
349 
350  bool cutsPassed = true;
351  if (scfg_.cut_h_apply_)
352  cutsPassed &= cv_h;
353  if (scfg_.cut_v_apply_)
354  cutsPassed &= cv_v;
355 
356  if (cutsPassed) {
357  pairsSelected++;
358 
359  h_q_cut_h_aft->Fill(cq_h);
360  h_q_cut_v_aft->Fill(cq_v);
361 
362  h2_cut_h_aft->Fill(trUp.x(), trDw.x());
363  h2_cut_v_aft->Fill(trUp.y(), trDw.y());
364 
365  m_h2_y_vs_x_aft_sel[scfg_.rp_N_.id_]->Fill(trUp.x(), trUp.y());
366  m_h2_y_vs_x_aft_sel[scfg_.rp_F_.id_]->Fill(trDw.x(), trDw.y());
367 
368  p_x_diffFN_vs_x_N->Fill(trUp.x(), trDw.x() - trUp.x());
369 
370  int idx = (trUp.x() - scfg_.rp_N_.x_slice_min_) / scfg_.rp_N_.x_slice_w_;
371  if (idx >= 0 && idx < scfg_.rp_N_.x_slice_n_) {
372  x_slice_plots_N[idx].fill(trUp.y(), trDw.y() - trUp.y(), debug);
373  }
374 
375  idx = (trDw.x() - scfg_.rp_F_.x_slice_min_) / scfg_.rp_F_.x_slice_w_;
376  if (idx >= 0 && idx < scfg_.rp_F_.x_slice_n_) {
377  x_slice_plots_F[idx].fill(trDw.y(), trDw.y() - trUp.y(), debug);
378  }
379  }
380  }
381  }
382 
383  return pairsSelected;
384 }
std::map< unsigned int, SlicePlots > x_slice_plots_F
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_bef_sel
uint32_t arm() const
Definition: CTPPSDetId.h:57
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_aft_sel
void Fill(long long x)
PPSAlignmentConfiguration::SectorConfig scfg_
#define debug
Definition: HDRShower.cc:19
auto const & tracks
cannot be loose
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
std::map< unsigned int, SlicePlots > x_slice_plots_N
std::map< unsigned int, MonitorElement * > m_h2_y_vs_x_mlt_sel

Member Data Documentation

◆ h2_cut_h_aft

MonitorElement* PPSAlignmentWorker::SectorData::h2_cut_h_aft

Definition at line 58 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ h2_cut_h_bef

MonitorElement* PPSAlignmentWorker::SectorData::h2_cut_h_bef

Definition at line 57 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ h2_cut_v_aft

MonitorElement* PPSAlignmentWorker::SectorData::h2_cut_v_aft

Definition at line 63 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ h2_cut_v_bef

MonitorElement* PPSAlignmentWorker::SectorData::h2_cut_v_bef

Definition at line 62 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ h_q_cut_h_aft

MonitorElement* PPSAlignmentWorker::SectorData::h_q_cut_h_aft

Definition at line 56 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ h_q_cut_h_bef

MonitorElement* PPSAlignmentWorker::SectorData::h_q_cut_h_bef

Definition at line 55 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ h_q_cut_v_aft

MonitorElement* PPSAlignmentWorker::SectorData::h_q_cut_v_aft

Definition at line 61 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ h_q_cut_v_bef

MonitorElement* PPSAlignmentWorker::SectorData::h_q_cut_v_bef

Definition at line 60 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ m_h2_y_vs_x_aft_sel

std::map<unsigned int, MonitorElement*> PPSAlignmentWorker::SectorData::m_h2_y_vs_x_aft_sel

Definition at line 52 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ m_h2_y_vs_x_bef_sel

std::map<unsigned int, MonitorElement*> PPSAlignmentWorker::SectorData::m_h2_y_vs_x_bef_sel

Definition at line 50 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ m_h2_y_vs_x_mlt_sel

std::map<unsigned int, MonitorElement*> PPSAlignmentWorker::SectorData::m_h2_y_vs_x_mlt_sel

Definition at line 51 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ p_x_diffFN_vs_x_N

MonitorElement* PPSAlignmentWorker::SectorData::p_x_diffFN_vs_x_N

Definition at line 66 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ p_y_diffFN_vs_y_F

MonitorElement* PPSAlignmentWorker::SectorData::p_y_diffFN_vs_y_F

Definition at line 67 of file PPSAlignmentWorker.cc.

◆ scfg_

PPSAlignmentConfiguration::SectorConfig PPSAlignmentWorker::SectorData::scfg_

Definition at line 47 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ x_slice_plots_F

std::map<unsigned int, SlicePlots> PPSAlignmentWorker::SectorData::x_slice_plots_F

Definition at line 79 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().

◆ x_slice_plots_N

std::map<unsigned int, SlicePlots> PPSAlignmentWorker::SectorData::x_slice_plots_N

Definition at line 79 of file PPSAlignmentWorker.cc.

Referenced by init(), and process().