CMS 3D CMS Logo

CTPPSProtonReconstructionValidator.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
11 
14 
16 
19 
23 
24 #include "TFile.h"
25 #include "TH1D.h"
26 
27 #include <map>
28 #include <string>
29 
30 //----------------------------------------------------------------------------------------------------
31 
33 public:
35 
36 private:
37  void analyze(const edm::Event&, const edm::EventSetup&) override;
38  void endJob() override;
39 
41  double chiSqCut_;
43 
44  struct RPPlots {
45  std::unique_ptr<TH1D> h_de_x, h_de_y;
47  : h_de_x(new TH1D("", ";#Deltax (mm)", 1000, -1., +1.)),
48  h_de_y(new TH1D("", ";#Deltay (mm)", 1000, -1., +1.)) {}
49 
50  void fill(double de_x, double de_y) {
51  h_de_x->Fill(de_x);
52  h_de_y->Fill(de_y);
53  }
54 
55  void write() const {
56  h_de_x->Write("h_de_x");
57  h_de_y->Write("h_de_y");
58  }
59  };
60 
61  std::map<unsigned int, RPPlots> rp_plots_;
62 };
63 
64 //----------------------------------------------------------------------------------------------------
65 
66 using namespace std;
67 using namespace edm;
68 
69 //----------------------------------------------------------------------------------------------------
70 
72  : tokenRecoProtons_(consumes<reco::ForwardProtonCollection>(iConfig.getParameter<edm::InputTag>("tagRecoProtons"))),
73  chiSqCut_(iConfig.getParameter<double>("chiSqCut")),
74  outputFile_(iConfig.getParameter<string>("outputFile")) {}
75 
76 //----------------------------------------------------------------------------------------------------
77 
79  // get conditions
81  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(hOpticalFunctions);
82 
83  // stop if conditions invalid
84  if (hOpticalFunctions->empty())
85  return;
86 
87  // get input
89  iEvent.getByToken(tokenRecoProtons_, hRecoProtons);
90 
91  // process tracks
92  for (const auto& pr : *hRecoProtons) {
93  if (!pr.validFit())
94  continue;
95 
96  if (pr.chi2() > chiSqCut_)
97  continue;
98 
99  for (const auto& tr : pr.contributingLocalTracks()) {
100  CTPPSDetId rpId(tr->rpId());
101  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
102 
103  // skip other than tracking RPs
104  if (rpId.subdetId() != CTPPSDetId::sdTrackingStrip && rpId.subdetId() != CTPPSDetId::sdTrackingPixel)
105  continue;
106 
107  // try to get optics for the RP
108  auto it = hOpticalFunctions->find(rpId);
109  if (it == hOpticalFunctions->end())
110  continue;
111 
112  // do propagation
113  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in_beam = {0., 0., 0., 0., 0.};
115  it->second.transport(k_in_beam, k_out_beam);
116 
118  -pr.vx(), -pr.thetaX(), pr.vy(), pr.thetaY(), pr.xi()}; // conversions: CMS --> LHC convention
120  it->second.transport(k_in, k_out);
121 
122  // fill plots
123  const double de_x = (k_out.x - k_out_beam.x) * 10. - tr->x(); // conversions: cm --> mm
124  const double de_y = (k_out.y - k_out_beam.y) * 10. - tr->y(); // conversions: cm --> mm
125 
126  rp_plots_[rpDecId].fill(de_x, de_y);
127  }
128  }
129 }
130 
131 //----------------------------------------------------------------------------------------------------
132 
134  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
135 
136  for (const auto& p : rp_plots_) {
137  gDirectory = f_out->mkdir(Form("%u", p.first));
138  p.second.write();
139  }
140 }
141 
142 //----------------------------------------------------------------------------------------------------
143 
EDAnalyzer.h
CTPPSProtonReconstructionValidator::endJob
void endJob() override
Definition: CTPPSProtonReconstructionValidator.cc:132
LHCInterpolatedOpticalFunctionsSet::Kinematics
proton kinematics description
Definition: LHCInterpolatedOpticalFunctionsSet.h:28
ESHandle.h
CTPPSProtonReconstructionValidator::RPPlots::fill
void fill(double de_x, double de_y)
Definition: CTPPSProtonReconstructionValidator.cc:51
edm::EDGetTokenT< reco::ForwardProtonCollection >
CTPPSProtonReconstructionValidator::RPPlots::RPPlots
RPPlots()
Definition: CTPPSProtonReconstructionValidator.cc:47
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
CTPPSLocalTrackLite.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
year_2016_postTS2_cff.rpId
rpId
Definition: year_2016_postTS2_cff.py:23
edm::Handle
Definition: AssociativeIterator.h:50
CTPPSDetId::sdTrackingStrip
Definition: CTPPSDetId.h:44
MakerMacros.h
CTPPSProtonReconstructionValidator::CTPPSProtonReconstructionValidator
CTPPSProtonReconstructionValidator(const edm::ParameterSet &)
Definition: CTPPSProtonReconstructionValidator.cc:70
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CTPPSProtonReconstructionValidator::tokenRecoProtons_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtons_
Definition: CTPPSProtonReconstructionValidator.cc:41
LHCInterpolatedOpticalFunctionsSet::Kinematics::y
double y
Definition: LHCInterpolatedOpticalFunctionsSet.h:31
CTPPSProtonReconstructionValidator::rp_plots_
std::map< unsigned int, RPPlots > rp_plots_
Definition: CTPPSProtonReconstructionValidator.cc:62
edm::ESHandle< LHCInterpolatedOpticalFunctionsSetCollection >
LHCInterpolatedOpticalFunctionsSetCollection.h
ForwardProtonFwd.h
LHCInterpolatedOpticalFunctionsSet::Kinematics::x
double x
Definition: LHCInterpolatedOpticalFunctionsSet.h:29
CTPPSDetId::sdTrackingPixel
Definition: CTPPSDetId.h:44
CTPPSProtonReconstructionValidator::RPPlots::write
void write() const
Definition: CTPPSProtonReconstructionValidator.cc:56
CTPPSProtonReconstructionValidator::RPPlots::h_de_x
std::unique_ptr< TH1D > h_de_x
Definition: CTPPSProtonReconstructionValidator.cc:46
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
CTPPSProtonReconstructionValidator
Definition: CTPPSProtonReconstructionValidator.cc:31
reco::ForwardProtonCollection
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
Definition: ForwardProtonFwd.h:25
edm::ParameterSet
Definition: ParameterSet.h:36
sipixeldigitoraw
Definition: SiPixelDigiToRaw.cc:38
Event.h
ForwardProton.h
CTPPSDetId
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:31
iEvent
int iEvent
Definition: GenABIO.cc:224
CTPPSInterpolatedOpticsRcd.h
CTPPSProtonReconstructionValidator::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: CTPPSProtonReconstructionValidator.cc:77
edm::EventSetup
Definition: EventSetup.h:57
alignCSCRings.de_y
de_y
Definition: alignCSCRings.py:86
get
#define get
CTPPSProtonReconstructionValidator::chiSqCut_
double chiSqCut_
Definition: CTPPSProtonReconstructionValidator.cc:42
CTPPSProtonReconstructionValidator::outputFile_
std::string outputFile_
Definition: CTPPSProtonReconstructionValidator.cc:43
std
Definition: JetResolutionObject.h:76
Frameworkfwd.h
ESWatcher.h
alignCSCRings.de_x
de_x
Definition: alignCSCRings.py:85
EventSetup.h
CTPPSDetId.h
CTPPSProtonReconstructionValidator::RPPlots
Definition: CTPPSProtonReconstructionValidator.cc:45
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
CTPPSInterpolatedOpticsRcd
Definition: CTPPSInterpolatedOpticsRcd.h:13
CTPPSProtonReconstructionValidator::RPPlots::h_de_y
std::unique_ptr< TH1D > h_de_y
Definition: CTPPSProtonReconstructionValidator.cc:46