CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
CTPPSProtonReconstructionPlotter Class Reference
Inheritance diagram for CTPPSProtonReconstructionPlotter:
edm::one::EDAnalyzer<> edm::one::EDAnalyzerBase edm::EDConsumerBase

Classes

struct  ArmCorrelationPlots
 
struct  AssociationCuts
 
struct  MultiRPPlots
 
struct  SingleMultiCorrelationPlots
 
struct  SingleRPPlots
 

Public Member Functions

 CTPPSProtonReconstructionPlotter (const edm::ParameterSet &)
 
 ~CTPPSProtonReconstructionPlotter () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer<>
 EDAnalyzer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void CalculateTimingTrackingDistance (const reco::ForwardProton &proton, const CTPPSLocalTrackLite &tr, const CTPPSGeometry &geometry, double &de_x, double &de_x_unc)
 
void endJob () override
 

Static Private Member Functions

static void profileToRMSGraph (TProfile *p, TGraphErrors *g)
 

Private Attributes

std::map< unsigned int, ArmCorrelationPlotsarmCorrelationPlots_
 
std::map< unsigned int, AssociationCutsassociation_cuts_
 
signed int maxNonEmptyEvents_
 
std::map< unsigned int, MultiRPPlotsmultiRPPlots_
 
signed int n_non_empty_events_
 
std::string outputFile_
 
std::unique_ptr< TProfile > p_x_L_diffNF_vs_x_L_N_
 
std::unique_ptr< TProfile > p_x_R_diffNF_vs_x_R_N_
 
std::unique_ptr< TProfile > p_y_L_diffNF_vs_y_L_N_
 
std::unique_ptr< TProfile > p_y_R_diffNF_vs_y_R_N_
 
unsigned int rpId_45_F_
 
unsigned int rpId_45_N_
 
unsigned int rpId_56_F_
 
unsigned int rpId_56_N_
 
std::map< unsigned int, SingleMultiCorrelationPlotssingleMultiCorrelationPlots_
 
std::map< unsigned int, SingleRPPlotssingleRPPlots_
 
edm::EDGetTokenT< reco::ForwardProtonCollectiontokenRecoProtonsMultiRP_
 
edm::EDGetTokenT< reco::ForwardProtonCollectiontokenRecoProtonsSingleRP_
 
edm::EDGetTokenT< CTPPSLocalTrackLiteCollectiontokenTracks_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 32 of file CTPPSProtonReconstructionPlotter.cc.

Constructor & Destructor Documentation

CTPPSProtonReconstructionPlotter::CTPPSProtonReconstructionPlotter ( const edm::ParameterSet ps)
explicit

Definition at line 462 of file CTPPSProtonReconstructionPlotter.cc.

References association_cuts_, edm::ParameterSet::getParameterSet(), and AlCaHLTBitMon_QueryRunRegistry::string.

463  : tokenTracks_(consumes<CTPPSLocalTrackLiteCollection>(ps.getParameter<edm::InputTag>("tagTracks"))),
465  consumes<reco::ForwardProtonCollection>(ps.getParameter<InputTag>("tagRecoProtonsSingleRP"))),
467  consumes<reco::ForwardProtonCollection>(ps.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
468 
469  rpId_45_N_(ps.getParameter<unsigned int>("rpId_45_N")),
470  rpId_45_F_(ps.getParameter<unsigned int>("rpId_45_F")),
471  rpId_56_N_(ps.getParameter<unsigned int>("rpId_56_N")),
472  rpId_56_F_(ps.getParameter<unsigned int>("rpId_56_F")),
473 
474  outputFile_(ps.getParameter<string>("outputFile")),
475  maxNonEmptyEvents_(ps.getUntrackedParameter<signed int>("maxNonEmptyEvents", -1)),
476 
477  p_x_L_diffNF_vs_x_L_N_(new TProfile("p_x_L_diffNF_vs_x_L_N", ";x_{LN};x_{LF} - x_{LN}", 100, 0., +20.)),
478  p_x_R_diffNF_vs_x_R_N_(new TProfile("p_x_R_diffNF_vs_x_R_N", ";x_{RN};x_{RF} - x_{RN}", 100, 0., +20.)),
479 
480  p_y_L_diffNF_vs_y_L_N_(new TProfile("p_y_L_diffNF_vs_y_L_N", ";y_{LN};y_{LF} - y_{LN}", 100, -20., +20.)),
481  p_y_R_diffNF_vs_y_R_N_(new TProfile("p_y_R_diffNF_vs_y_R_N", ";y_{RN};y_{RF} - y_{RN}", 100, -20., +20.)),
482 
484  for (const std::string &sector : {"45", "56"}) {
485  const unsigned int arm = (sector == "45") ? 0 : 1;
486  association_cuts_[arm].load(ps.getParameterSet("association_cuts_" + sector));
487  }
488 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tokenTracks_
ParameterSet const & getParameterSet(std::string const &) const
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
std::map< unsigned int, AssociationCuts > association_cuts_
CTPPSProtonReconstructionPlotter::~CTPPSProtonReconstructionPlotter ( )
inlineoverride

Definition at line 35 of file CTPPSProtonReconstructionPlotter.cc.

References analyze(), and endJob().

35 {}

Member Function Documentation

void CTPPSProtonReconstructionPlotter::analyze ( const edm::Event event,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 519 of file CTPPSProtonReconstructionPlotter.cc.

References CTPPSDetId::arm(), armCorrelationPlots_, association_cuts_, CalculateTimingTrackingDistance(), reco::ForwardProton::contributingLocalTracks(), alignCSCRings::de_x, Exception, edm::EventSetup::get(), CTPPSLocalTrackLite::getX(), CTPPSLocalTrackLite::getY(), training_settings::idx, match(), maxNonEmptyEvents_, multiRPPlots_, n_non_empty_events_, p_x_L_diffNF_vs_x_L_N_, p_x_R_diffNF_vs_x_R_N_, p_y_L_diffNF_vs_y_L_N_, p_y_R_diffNF_vs_y_R_N_, year_2016_postTS2_cff::rpId, rpId_45_F_, rpId_45_N_, rpId_56_F_, rpId_56_N_, CTPPSDetId::sdTrackingPixel, CTPPSDetId::sdTrackingStrip, singleMultiCorrelationPlots_, singleRPPlots_, tokenRecoProtonsMultiRP_, tokenRecoProtonsSingleRP_, and tokenTracks_.

Referenced by ~CTPPSProtonReconstructionPlotter().

519  {
520  // get input
522  event.getByToken(tokenTracks_, hTracks);
523 
524  Handle<reco::ForwardProtonCollection> hRecoProtonsSingleRP;
525  event.getByToken(tokenRecoProtonsSingleRP_, hRecoProtonsSingleRP);
526 
527  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
528  event.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
529 
530  if (!hRecoProtonsSingleRP->empty())
532 
534  throw cms::Exception("CTPPSProtonReconstructionPlotter") << "Number of non empty events reached maximum.";
535 
536  // get conditions
538  iSetup.get<VeryForwardRealGeometryRecord>().get(hGeometry);
539 
540  // track plots
541  const CTPPSLocalTrackLite *tr_L_N = nullptr;
542  const CTPPSLocalTrackLite *tr_L_F = nullptr;
543  const CTPPSLocalTrackLite *tr_R_N = nullptr;
544  const CTPPSLocalTrackLite *tr_R_F = nullptr;
545  std::map<unsigned int, unsigned int> armTrackCounter, armTimingTrackCounter;
546  std::map<unsigned int, unsigned int> armTrackCounter_N, armTrackCounter_F;
547 
548  for (const auto &tr : *hTracks) {
549  CTPPSDetId rpId(tr.getRPId());
550  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
551 
552  if (decRPId == rpId_45_N_) {
553  tr_L_N = &tr;
554  armTrackCounter_N[0]++;
555  }
556  if (decRPId == rpId_45_F_) {
557  tr_L_F = &tr;
558  armTrackCounter_F[0]++;
559  }
560  if (decRPId == rpId_56_N_) {
561  tr_R_N = &tr;
562  armTrackCounter_N[1]++;
563  }
564  if (decRPId == rpId_56_F_) {
565  tr_R_F = &tr;
566  armTrackCounter_F[1]++;
567  }
568 
569  armTrackCounter[rpId.arm()]++;
570 
571  const bool trackerRP =
572  (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
573  if (!trackerRP)
574  armTimingTrackCounter[rpId.arm()]++;
575  }
576 
577  if (tr_L_N && tr_L_F) {
578  p_x_L_diffNF_vs_x_L_N_->Fill(tr_L_N->getX(), tr_L_F->getX() - tr_L_N->getX());
579  p_y_L_diffNF_vs_y_L_N_->Fill(tr_L_N->getY(), tr_L_F->getY() - tr_L_N->getY());
580  }
581 
582  if (tr_R_N && tr_R_F) {
583  p_x_R_diffNF_vs_x_R_N_->Fill(tr_R_N->getX(), tr_R_F->getX() - tr_R_N->getX());
584  p_y_R_diffNF_vs_y_R_N_->Fill(tr_R_N->getY(), tr_R_F->getY() - tr_R_N->getY());
585  }
586 
587  // initialise multiplicity counters
588  std::map<unsigned int, unsigned int> singleRPMultiplicity, multiRPMultiplicity;
589  singleRPMultiplicity[rpId_45_N_] = singleRPMultiplicity[rpId_45_F_] = singleRPMultiplicity[rpId_56_N_] =
590  singleRPMultiplicity[rpId_56_F_] = 0;
591  multiRPMultiplicity[0] = multiRPMultiplicity[1] = 0;
592 
593  // make single-RP-reco plots
594  for (const auto &proton : *hRecoProtonsSingleRP) {
595  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
596  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
597 
598  const bool n1f1 = (armTrackCounter_N[rpId.arm()] == 1 && armTrackCounter_F[rpId.arm()] == 1);
599 
600  singleRPPlots_[decRPId].fill(proton, armTrackCounter[rpId.arm()], n1f1);
601 
602  if (proton.validFit())
603  singleRPMultiplicity[decRPId]++;
604  }
605 
606  for (const auto it : singleRPMultiplicity)
607  singleRPPlots_[it.first].h_multiplicity->Fill(it.second);
608 
609  // make multi-RP-reco plots
610  for (const auto &proton : *hRecoProtonsMultiRP) {
611  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
612  unsigned int armId = rpId.arm();
613 
614  const bool n1f1 = (armTrackCounter_N[armId] == 1 && armTrackCounter_F[armId] == 1);
615 
616  multiRPPlots_[armId].fill(proton, armTrackCounter[armId], n1f1);
617 
618  if (proton.validFit())
619  multiRPMultiplicity[armId]++;
620  }
621 
622  for (const auto it : multiRPMultiplicity) {
623  const auto &pl = multiRPPlots_[it.first];
624  pl.h_multiplicity->Fill(it.second);
625  pl.h2_timing_tracks_vs_prot_mult->Fill(it.second, armTimingTrackCounter[it.first]);
626  }
627 
628  // define "clean condition" for each arm
629  map<unsigned int, bool> clCo;
630  clCo[0] = (singleRPMultiplicity[rpId_45_N_] && singleRPMultiplicity[rpId_45_F_] && multiRPMultiplicity[0] == 1);
631  clCo[1] = (singleRPMultiplicity[rpId_56_N_] && singleRPMultiplicity[rpId_56_F_] && multiRPMultiplicity[1] == 1);
632 
633  // plot distances between multi-RP protons and timing tracks in the same arm
634  for (const auto &proton : *hRecoProtonsMultiRP) {
635  if (!proton.validFit())
636  continue;
637 
638  CTPPSDetId rpId_proton((*proton.contributingLocalTracks().begin())->getRPId());
639  unsigned int armId = rpId_proton.arm();
640  const auto &pl = multiRPPlots_[armId];
641 
642  for (const auto &tr : *hTracks) {
643  CTPPSDetId rpId_tr(tr.getRPId());
644  if (rpId_tr.arm() != armId)
645  continue;
646 
647  const bool trackerRP =
648  (rpId_tr.subdetId() == CTPPSDetId::sdTrackingStrip || rpId_tr.subdetId() == CTPPSDetId::sdTrackingPixel);
649  if (trackerRP)
650  continue;
651 
652  double de_x = 0., de_x_unc = 0.;
653  CalculateTimingTrackingDistance(proton, tr, *hGeometry, de_x, de_x_unc);
654 
655  const double rd = (de_x_unc > 0.) ? de_x / de_x_unc : -1E10;
656  const auto &ac = association_cuts_[armId];
657  const bool match = (ac.ti_tr_min <= fabs(rd) && fabs(rd) <= ac.ti_tr_max);
658 
659  pl.h_de_x_timing_vs_tracking->Fill(de_x);
660  pl.h_de_x_rel_timing_vs_tracking->Fill(rd);
661  pl.h_de_x_match_timing_vs_tracking->Fill(match ? 1. : 0.);
662 
663  if (clCo[armId] && armTimingTrackCounter[armId] == 1) {
664  pl.h_de_x_timing_vs_tracking_ClCo->Fill(de_x);
665  pl.h_de_x_rel_timing_vs_tracking_ClCo->Fill(rd);
666  pl.h_de_x_match_timing_vs_tracking_ClCo->Fill(match ? 1. : 0.);
667  }
668  }
669  }
670 
671  // plot xy maps
672  for (const auto &proton : *hRecoProtonsMultiRP) {
673  if (!proton.validFit())
674  continue;
675 
676  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
677  unsigned int armId = rpId.arm();
678  const auto &pl = multiRPPlots_[armId];
679  const auto &nTimingTracks = armTimingTrackCounter[armId];
680 
681  if (!clCo[armId])
682  continue;
683 
684  double x_ref = 0., y_ref = 0.;
685  if (armId == 0) {
686  x_ref = tr_L_N->getX();
687  y_ref = tr_L_N->getY();
688  }
689  if (armId == 1) {
690  x_ref = tr_R_N->getX();
691  y_ref = tr_R_N->getY();
692  }
693 
694  if (nTimingTracks == 0)
695  pl.h2_y_vs_x_tt0_ClCo->Fill(x_ref, y_ref);
696  if (nTimingTracks == 1)
697  pl.h2_y_vs_x_tt1_ClCo->Fill(x_ref, y_ref);
698  if (nTimingTracks > 1)
699  pl.h2_y_vs_x_ttm_ClCo->Fill(x_ref, y_ref);
700  }
701 
702  // make correlation plots
703  for (const auto &proton_m : *hRecoProtonsMultiRP) {
704  CTPPSDetId rpId_m((*proton_m.contributingLocalTracks().begin())->getRPId());
705  unsigned int arm = rpId_m.arm();
706 
707  const reco::ForwardProton *p_s_N = nullptr, *p_s_F = nullptr;
708 
709  for (const auto &proton_s : *hRecoProtonsSingleRP) {
710  // skip if source of single-RP reco not included in multi-RP reco
711  const auto key_s = proton_s.contributingLocalTracks()[0].key();
712  bool compatible = false;
713  for (const auto &tr_m : proton_m.contributingLocalTracks()) {
714  if (tr_m.key() == key_s) {
715  compatible = true;
716  break;
717  }
718  }
719 
720  if (!compatible)
721  continue;
722 
723  // fill single-multi plots
724  CTPPSDetId rpId_s((*proton_s.contributingLocalTracks().begin())->getRPId());
725  const unsigned int idx = rpId_s.arm() * 1000 + rpId_s.station() * 100 + rpId_s.rp() * 10 + rpId_s.arm();
726  singleMultiCorrelationPlots_[idx].fill(proton_s, proton_m);
727 
728  // select protons for arm-correlation plots
729  const unsigned int rpDecId_s = rpId_s.arm() * 100 + rpId_s.station() * 10 + rpId_s.rp();
730  if (rpDecId_s == rpId_45_N_ || rpDecId_s == rpId_56_N_)
731  p_s_N = &proton_s;
732  if (rpDecId_s == rpId_45_F_ || rpDecId_s == rpId_56_F_)
733  p_s_F = &proton_s;
734  }
735 
736  // fill arm-correlation plots
737  if (p_s_N && p_s_F)
738  armCorrelationPlots_[arm].fill(*p_s_N, *p_s_F, proton_m);
739  }
740 }
std::map< unsigned int, MultiRPPlots > multiRPPlots_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
std::map< unsigned int, SingleRPPlots > singleRPPlots_
Local (=single RP) track with essential information only.
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tokenTracks_
float getX() const
returns the horizontal track position
Event setup record containing the real (actual) geometry information.
float getY() const
returns the vertical track position
uint32_t arm() const
Definition: CTPPSDetId.h:51
const CTPPSLocalTrackLiteRefVector & contributingLocalTracks() const
list of RP tracks that contributed to this global track
std::map< unsigned int, ArmCorrelationPlots > armCorrelationPlots_
void CalculateTimingTrackingDistance(const reco::ForwardProton &proton, const CTPPSLocalTrackLite &tr, const CTPPSGeometry &geometry, double &de_x, double &de_x_unc)
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
std::map< unsigned int, SingleMultiCorrelationPlots > singleMultiCorrelationPlots_
T get() const
Definition: EventSetup.h:71
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
std::map< unsigned int, AssociationCuts > association_cuts_
void CTPPSProtonReconstructionPlotter::CalculateTimingTrackingDistance ( const reco::ForwardProton proton,
const CTPPSLocalTrackLite tr,
const CTPPSGeometry geometry,
double &  de_x,
double &  de_x_unc 
)
private

Definition at line 492 of file CTPPSProtonReconstructionPlotter.cc.

References edm::RefVector< C, T, F >::at(), reco::ForwardProton::contributingLocalTracks(), CTPPSLocalTrackLite::getRPId(), CTPPSGeometry::getRPTranslation(), CTPPSLocalTrackLite::getX(), CTPPSLocalTrackLite::getXUnc(), mathSSE::sqrt(), and z.

Referenced by analyze(), and profileToRMSGraph().

496  {
497  // identify tracking-RP tracks
498  const auto &tr_i = *proton.contributingLocalTracks().at(0);
499  const auto &tr_j = *proton.contributingLocalTracks().at(1);
500 
501  const double z_i = geometry.getRPTranslation(tr_i.getRPId()).z();
502  const double z_j = geometry.getRPTranslation(tr_j.getRPId()).z();
503 
504  // interpolation from tracking RPs
505  const double z_ti =
506  -geometry.getRPTranslation(tr_ti.getRPId()).z(); // the minus sign fixes a bug in the diamond geometry
507  const double f_i = (z_ti - z_j) / (z_i - z_j), f_j = (z_i - z_ti) / (z_i - z_j);
508  const double x_inter = f_i * tr_i.getX() + f_j * tr_j.getX();
509  const double x_inter_unc_sq =
510  f_i * f_i * tr_i.getXUnc() * tr_i.getXUnc() + f_j * f_j * tr_j.getXUnc() * tr_j.getXUnc();
511 
512  // save distance
513  de_x = tr_ti.getX() - x_inter;
514  de_x_unc = sqrt(tr_ti.getXUnc() * tr_ti.getXUnc() + x_inter_unc_sq);
515 }
T sqrt(T t)
Definition: SSEVec.h:18
const CTPPSLocalTrackLiteRefVector & contributingLocalTracks() const
list of RP tracks that contributed to this global track
CLHEP::Hep3Vector getRPTranslation(unsigned int id) const
value_type const at(size_type idx) const
Retrieve an element of the RefVector.
Definition: RefVector.h:88
void CTPPSProtonReconstructionPlotter::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 744 of file CTPPSProtonReconstructionPlotter.cc.

References armCorrelationPlots_, DEFINE_FWK_MODULE, multiRPPlots_, n_non_empty_events_, outputFile_, p_x_L_diffNF_vs_x_L_N_, p_x_R_diffNF_vs_x_R_N_, p_y_L_diffNF_vs_y_L_N_, p_y_R_diffNF_vs_y_R_N_, singleMultiCorrelationPlots_, and singleRPPlots_.

Referenced by ~CTPPSProtonReconstructionPlotter().

744  {
745  LogInfo("CTPPSProtonReconstructionPlotter") << "n_non_empty_events = " << n_non_empty_events_;
746 
747  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
748 
749  p_x_L_diffNF_vs_x_L_N_->Write();
750  p_x_R_diffNF_vs_x_R_N_->Write();
751 
752  p_y_L_diffNF_vs_y_L_N_->Write();
753  p_y_R_diffNF_vs_y_R_N_->Write();
754 
755  TDirectory *d_singleRPPlots = f_out->mkdir("singleRPPlots");
756  for (const auto &it : singleRPPlots_) {
757  gDirectory = d_singleRPPlots->mkdir(Form("rp%u", it.first));
758  it.second.write();
759  }
760 
761  TDirectory *d_multiRPPlots = f_out->mkdir("multiRPPlots");
762  for (const auto &it : multiRPPlots_) {
763  gDirectory = d_multiRPPlots->mkdir(Form("arm%u", it.first));
764  it.second.write();
765  }
766 
767  TDirectory *d_singleMultiCorrelationPlots = f_out->mkdir("singleMultiCorrelationPlots");
768  for (const auto &it : singleMultiCorrelationPlots_) {
769  unsigned int si_rp = it.first / 10;
770  unsigned int mu_arm = it.first % 10;
771 
772  gDirectory = d_singleMultiCorrelationPlots->mkdir(Form("si_rp%u_mu_arm%u", si_rp, mu_arm));
773  it.second.write();
774  }
775 
776  TDirectory *d_armCorrelationPlots = f_out->mkdir("armCorrelationPlots");
777  for (const auto &it : armCorrelationPlots_) {
778  gDirectory = d_armCorrelationPlots->mkdir(Form("arm%u", it.first));
779  it.second.write();
780  }
781 }
std::map< unsigned int, MultiRPPlots > multiRPPlots_
std::map< unsigned int, SingleRPPlots > singleRPPlots_
std::map< unsigned int, ArmCorrelationPlots > armCorrelationPlots_
std::map< unsigned int, SingleMultiCorrelationPlots > singleMultiCorrelationPlots_
static void CTPPSProtonReconstructionPlotter::profileToRMSGraph ( TProfile *  p,
TGraphErrors *  g 
)
inlinestaticprivate

Definition at line 65 of file CTPPSProtonReconstructionPlotter.cc.

References EnergyCorrector::c, CalculateTimingTrackingDistance(), alignCSCRings::de_x, training_settings::idx, N, and mathSSE::sqrt().

Referenced by CTPPSProtonReconstructionPlotter::MultiRPPlots::write().

65  {
66  for (int bi = 1; bi <= p->GetNbinsX(); ++bi) {
67  double c = p->GetBinCenter(bi);
68 
69  double N = p->GetBinEntries(bi);
70  double Sy = p->GetBinContent(bi) * N;
71  double Syy = p->GetSumw2()->At(bi);
72 
73  double si_sq = Syy / N - Sy * Sy / N / N;
74  double si = (si_sq >= 0.) ? sqrt(si_sq) : 0.;
75  double si_unc_sq = si_sq / 2. / N; // Gaussian approximation
76  double si_unc = (si_unc_sq >= 0.) ? sqrt(si_unc_sq) : 0.;
77 
78  int idx = g->GetN();
79  g->SetPoint(idx, c, si);
80  g->SetPointError(idx, 0., si_unc);
81  }
82  }
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
T sqrt(T t)
Definition: SSEVec.h:18
#define N
Definition: blowfish.cc:9

Member Data Documentation

std::map<unsigned int, ArmCorrelationPlots> CTPPSProtonReconstructionPlotter::armCorrelationPlots_
private

Definition at line 447 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

std::map<unsigned int, AssociationCuts> CTPPSProtonReconstructionPlotter::association_cuts_
private
signed int CTPPSProtonReconstructionPlotter::maxNonEmptyEvents_
private

Definition at line 63 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().

std::map<unsigned int, MultiRPPlots> CTPPSProtonReconstructionPlotter::multiRPPlots_
private

Definition at line 356 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

signed int CTPPSProtonReconstructionPlotter::n_non_empty_events_
private

Definition at line 452 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

std::string CTPPSProtonReconstructionPlotter::outputFile_
private

Definition at line 61 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by endJob().

std::unique_ptr<TProfile> CTPPSProtonReconstructionPlotter::p_x_L_diffNF_vs_x_L_N_
private

Definition at line 449 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

std::unique_ptr<TProfile> CTPPSProtonReconstructionPlotter::p_x_R_diffNF_vs_x_R_N_
private

Definition at line 449 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

std::unique_ptr<TProfile> CTPPSProtonReconstructionPlotter::p_y_L_diffNF_vs_y_L_N_
private

Definition at line 450 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

std::unique_ptr<TProfile> CTPPSProtonReconstructionPlotter::p_y_R_diffNF_vs_y_R_N_
private

Definition at line 450 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

unsigned int CTPPSProtonReconstructionPlotter::rpId_45_F_
private

Definition at line 46 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().

unsigned int CTPPSProtonReconstructionPlotter::rpId_45_N_
private

Definition at line 46 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().

unsigned int CTPPSProtonReconstructionPlotter::rpId_56_F_
private

Definition at line 47 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().

unsigned int CTPPSProtonReconstructionPlotter::rpId_56_N_
private

Definition at line 47 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().

std::map<unsigned int, SingleMultiCorrelationPlots> CTPPSProtonReconstructionPlotter::singleMultiCorrelationPlots_
private

Definition at line 405 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

std::map<unsigned int, SingleRPPlots> CTPPSProtonReconstructionPlotter::singleRPPlots_
private

Definition at line 148 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze(), and endJob().

edm::EDGetTokenT<reco::ForwardProtonCollection> CTPPSProtonReconstructionPlotter::tokenRecoProtonsMultiRP_
private

Definition at line 44 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().

edm::EDGetTokenT<reco::ForwardProtonCollection> CTPPSProtonReconstructionPlotter::tokenRecoProtonsSingleRP_
private

Definition at line 43 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().

edm::EDGetTokenT<CTPPSLocalTrackLiteCollection> CTPPSProtonReconstructionPlotter::tokenTracks_
private

Definition at line 42 of file CTPPSProtonReconstructionPlotter.cc.

Referenced by analyze().