CMS 3D CMS Logo

CTPPSProtonReconstructionDiffPlotter.cc
Go to the documentation of this file.
1 /****************************************************************************
2 * Authors:
3 * Jan Kašpar (jan.kaspar@gmail.com)
4 ****************************************************************************/
5 
6 // TODO: clean
15 
17 
21 
24 
25 #include "TFile.h"
26 #include "TGraphErrors.h"
27 #include "TH1D.h"
28 #include "TH2D.h"
29 #include "TProfile.h"
30 
31 //----------------------------------------------------------------------------------------------------
32 
34 public:
37 
38 private:
39  void analyze(const edm::Event &, const edm::EventSetup &) override;
40 
41  void endJob() override;
42 
45 
47 
48  std::unique_ptr<TH1D> h_de_xi, h_de_th_x, h_de_th_y, h_de_vtx_y;
49 };
50 
51 //----------------------------------------------------------------------------------------------------
52 
53 using namespace std;
54 using namespace edm;
55 
56 //----------------------------------------------------------------------------------------------------
57 
59  : tokenRecoProtonsRef_(consumes<reco::ForwardProtonCollection>(ps.getParameter<InputTag>("tagRecoProtonsRef"))),
60  tokenRecoProtonsTest_(consumes<reco::ForwardProtonCollection>(ps.getParameter<InputTag>("tagRecoProtonsTest"))),
61 
62  outputFile_(ps.getParameter<string>("outputFile")),
63 
64  h_de_xi(new TH1D("", ";#Delta#xi", 200, -0.01, +0.01)),
65  h_de_th_x(new TH1D("", ";#Delta#theta_{x}", 200, -100E-6, +100E-6)),
66  h_de_th_y(new TH1D("", ";#Delta#theta_{y}", 200, -100E-6, +100E-6)),
67  h_de_vtx_y(new TH1D("", ";#Deltay^{*} (cm)", 200, -0.01, +0.01)) {}
68 
69 //----------------------------------------------------------------------------------------------------
70 
72  // get input
74  event.getByToken(tokenRecoProtonsRef_, hRecoProtonsRef);
75 
77  event.getByToken(tokenRecoProtonsTest_, hRecoProtonsTest);
78 
79  if (hRecoProtonsRef->size() != hRecoProtonsTest->size()) {
80  edm::LogWarning("CTPPSProtonReconstructionDiffPlotter::analyze")
81  << "Different number of Ref and Test protons. Skipping event.";
82  return;
83  }
84 
85  for (unsigned int i = 0; i < hRecoProtonsRef->size(); ++i) {
86  const auto &pr_ref = hRecoProtonsRef->at(i);
87  const auto &pr_test = hRecoProtonsTest->at(i);
88 
89  if (!pr_ref.validFit() || !pr_test.validFit())
90  continue;
91 
92  h_de_xi->Fill(pr_test.xi() - pr_ref.xi());
93  h_de_th_x->Fill(pr_test.thetaX() - pr_ref.thetaX());
94  h_de_th_y->Fill(pr_test.thetaY() - pr_ref.thetaY());
95  h_de_vtx_y->Fill(pr_test.vy() - pr_ref.vy());
96  }
97 }
98 
99 //----------------------------------------------------------------------------------------------------
100 
102  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
103 
104  h_de_xi->Write("h_de_xi");
105  h_de_th_x->Write("h_de_th_x");
106  h_de_th_y->Write("h_de_th_y");
107  h_de_vtx_y->Write("h_de_vtx_y");
108 }
109 
110 //----------------------------------------------------------------------------------------------------
111 
EDAnalyzer.h
mps_fire.i
i
Definition: mps_fire.py:428
ESHandle.h
edm::EDGetTokenT< reco::ForwardProtonCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
CTPPSProtonReconstructionDiffPlotter::h_de_th_x
std::unique_ptr< TH1D > h_de_th_x
Definition: CTPPSProtonReconstructionDiffPlotter.cc:49
CTPPSProtonReconstructionDiffPlotter::outputFile_
std::string outputFile_
Definition: CTPPSProtonReconstructionDiffPlotter.cc:47
CTPPSLocalTrackLite.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
CTPPSProtonReconstructionDiffPlotter::~CTPPSProtonReconstructionDiffPlotter
~CTPPSProtonReconstructionDiffPlotter() override
Definition: CTPPSProtonReconstructionDiffPlotter.cc:37
edm::Handle
Definition: AssociativeIterator.h:50
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
CTPPSProtonReconstructionDiffPlotter::tokenRecoProtonsRef_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsRef_
Definition: CTPPSProtonReconstructionDiffPlotter.cc:44
CTPPSProtonReconstructionDiffPlotter::h_de_xi
std::unique_ptr< TH1D > h_de_xi
Definition: CTPPSProtonReconstructionDiffPlotter.cc:49
CTPPSGeometry.h
MakerMacros.h
CTPPSProtonReconstructionDiffPlotter::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: CTPPSProtonReconstructionDiffPlotter.cc:70
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ForwardProtonFwd.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
reco::ForwardProtonCollection
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
Definition: ForwardProtonFwd.h:25
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
CTPPSProtonReconstructionDiffPlotter::tokenRecoProtonsTest_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsTest_
Definition: CTPPSProtonReconstructionDiffPlotter.cc:45
ForwardProton.h
CTPPSProtonReconstructionDiffPlotter
Definition: CTPPSProtonReconstructionDiffPlotter.cc:32
edm::EventSetup
Definition: EventSetup.h:57
VeryForwardRealGeometryRecord.h
CTPPSProtonReconstructionDiffPlotter::h_de_vtx_y
std::unique_ptr< TH1D > h_de_vtx_y
Definition: CTPPSProtonReconstructionDiffPlotter.cc:49
std
Definition: JetResolutionObject.h:76
Frameworkfwd.h
CTPPSProtonReconstructionDiffPlotter::endJob
void endJob() override
Definition: CTPPSProtonReconstructionDiffPlotter.cc:100
CTPPSProtonReconstructionDiffPlotter::CTPPSProtonReconstructionDiffPlotter
CTPPSProtonReconstructionDiffPlotter(const edm::ParameterSet &)
Definition: CTPPSProtonReconstructionDiffPlotter.cc:57
EventSetup.h
CTPPSDetId.h
Exception.h
ParameterSet.h
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
CTPPSProtonReconstructionDiffPlotter::h_de_th_y
std::unique_ptr< TH1D > h_de_th_y
Definition: CTPPSProtonReconstructionDiffPlotter.cc:49