CMS 3D CMS Logo

CTPPSOpticsPlotter.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
10 
13 
15 
18 
19 #include "TFile.h"
20 #include "TGraph.h"
21 
22 #include <map>
23 #include <string>
24 
25 //----------------------------------------------------------------------------------------------------
26 
28 public:
29  explicit CTPPSOpticsPlotter(const edm::ParameterSet&);
30 
31 private:
32  void analyze(const edm::Event&, const edm::EventSetup&) override;
33  void endJob() override;
34 
36 
38 
39  struct RPPlots {
40  std::unique_ptr<TGraph> g_v_x_vs_xi, g_L_x_vs_xi, g_x_D_vs_xi;
41  std::unique_ptr<TGraph> g_v_y_vs_xi, g_L_y_vs_xi, g_y_D_vs_xi;
42  std::unique_ptr<TGraph> h_y_vs_x_disp;
43 
45  : g_v_x_vs_xi(new TGraph),
46  g_L_x_vs_xi(new TGraph),
47  g_x_D_vs_xi(new TGraph),
48  g_v_y_vs_xi(new TGraph),
49  g_L_y_vs_xi(new TGraph),
50  g_y_D_vs_xi(new TGraph),
51  h_y_vs_x_disp(new TGraph) {}
52 
53  void write() const {
54  g_v_x_vs_xi->SetTitle(";xi;v_{x}");
55  g_v_x_vs_xi->Write("g_v_x_vs_xi");
56 
57  g_L_x_vs_xi->SetTitle(";xi;L_{x} (cm)");
58  g_L_x_vs_xi->Write("g_L_x_vs_xi");
59 
60  g_x_D_vs_xi->SetTitle(";xi;x_{D} (cm)");
61  g_x_D_vs_xi->Write("g_x_D_vs_xi");
62 
63  g_v_y_vs_xi->SetTitle(";xi;v_{y}");
64  g_v_y_vs_xi->Write("g_v_y_vs_xi");
65 
66  g_L_y_vs_xi->SetTitle(";xi;L_{y} (cm)");
67  g_L_y_vs_xi->Write("g_L_y_vs_xi");
68 
69  g_y_D_vs_xi->SetTitle(";xi;y_{D} (cm)");
70  g_y_D_vs_xi->Write("g_y_D_vs_xi");
71 
72  h_y_vs_x_disp->SetTitle(";x (cm);y (cm)");
73  h_y_vs_x_disp->Write("h_y_vs_x_disp");
74  }
75  };
76 
77  std::map<unsigned int, RPPlots> rp_plots_;
78 };
79 
80 //----------------------------------------------------------------------------------------------------
81 
82 using namespace std;
83 using namespace edm;
84 
85 //----------------------------------------------------------------------------------------------------
86 
88  : opticsLabel_(iConfig.getParameter<std::string>("opticsLabel")),
89  outputFile_(iConfig.getParameter<string>("outputFile")) {}
90 
91 //----------------------------------------------------------------------------------------------------
92 
94  // stop if plots already made
95  if (!rp_plots_.empty())
96  return;
97 
98  // get conditions
100  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(opticsLabel_, hOpticalFunctions);
101 
102  // stop if conditions invalid
103  if (hOpticalFunctions->empty())
104  return;
105 
106  // make plots
107  for (const auto& it : *hOpticalFunctions) {
108  CTPPSDetId rpId(it.first);
109  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
110 
111  auto& pl = rp_plots_[rpDecId];
112 
113  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in_beam = {0., 0., 0., 0., 0.};
115  it.second.transport(k_in_beam, k_out_beam);
116 
117  const double vtx_ep = 1E-4; // cm
118  const double th_ep = 1E-6; // rad
119 
120  for (double xi = 0.; xi < 0.30001; xi += 0.001) {
121  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in_xi = {0., 0., 0., 0., xi};
123  it.second.transport(k_in_xi, k_out_xi);
124 
125  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in_xi_vtx_x = {vtx_ep, 0., 0., 0., xi};
127  it.second.transport(k_in_xi_vtx_x, k_out_xi_vtx_x);
128 
129  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in_xi_th_x = {0., th_ep, 0., 0., xi};
131  it.second.transport(k_in_xi_th_x, k_out_xi_th_x);
132 
133  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in_xi_vtx_y = {0., 0., vtx_ep, 0., xi};
135  it.second.transport(k_in_xi_vtx_y, k_out_xi_vtx_y);
136 
137  LHCInterpolatedOpticalFunctionsSet::Kinematics k_in_xi_th_y = {0., 0., 0., th_ep, xi};
139  it.second.transport(k_in_xi_th_y, k_out_xi_th_y);
140 
141  int idx = pl.g_v_x_vs_xi->GetN();
142 
143  pl.g_v_x_vs_xi->SetPoint(idx, xi, (k_out_xi_vtx_x.x - k_out_xi.x) / vtx_ep);
144  pl.g_L_x_vs_xi->SetPoint(idx, xi, (k_out_xi_th_x.x - k_out_xi.x) / th_ep);
145  pl.g_x_D_vs_xi->SetPoint(idx, xi, k_out_xi.x - k_out_beam.x);
146 
147  pl.g_v_y_vs_xi->SetPoint(idx, xi, (k_out_xi_vtx_y.y - k_out_xi.y) / vtx_ep);
148  pl.g_L_y_vs_xi->SetPoint(idx, xi, (k_out_xi_th_y.y - k_out_xi.y) / th_ep);
149  pl.g_y_D_vs_xi->SetPoint(idx, xi, k_out_xi.y - k_out_beam.y);
150 
151  pl.h_y_vs_x_disp->SetPoint(idx, k_out_xi.x - k_out_beam.x, k_out_xi.y - k_out_beam.y);
152  }
153  }
154 }
155 
156 //----------------------------------------------------------------------------------------------------
157 
159  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
160 
161  for (const auto& p : rp_plots_) {
162  gDirectory = f_out->mkdir(Form("%u", p.first));
163  p.second.write();
164  }
165 }
166 
167 //----------------------------------------------------------------------------------------------------
168 
uint32_t station() const
Definition: CTPPSDetId.h:58
void endJob() override
void analyze(const edm::Event &, const edm::EventSetup &) override
std::unique_ptr< TGraph > g_L_x_vs_xi
CTPPSOpticsPlotter(const edm::ParameterSet &)
std::unique_ptr< TGraph > g_x_D_vs_xi
std::unique_ptr< TGraph > g_L_y_vs_xi
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
uint32_t arm() const
Definition: CTPPSDetId.h:51
std::unique_ptr< TGraph > g_v_y_vs_xi
std::map< unsigned int, RPPlots > rp_plots_
std::unique_ptr< TGraph > g_v_x_vs_xi
std::unique_ptr< TGraph > h_y_vs_x_disp
std::unique_ptr< TGraph > g_y_D_vs_xi
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
HLT enums.
T get() const
Definition: EventSetup.h:71
uint32_t rp() const
Definition: CTPPSDetId.h:65