CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BSvsPVHistogramMaker.cc
Go to the documentation of this file.
9 #include "TH1F.h"
10 #include "TH2F.h"
11 #include "TProfile.h"
12 
14  : _currdir(nullptr),
15  m_maxLS(100),
16  useSlope_(true),
17  _trueOnly(true),
18  _runHisto(true),
19  _runHistoProfile(true),
20  _runHistoBXProfile(true),
21  _runHistoBX2D(false),
22  _histoParameters(),
23  _rhm(iC) {}
24 
26  : _currdir(nullptr),
27  m_maxLS(iConfig.getParameter<unsigned int>("maxLSBeforeRebin")),
28  useSlope_(iConfig.getParameter<bool>("useSlope")),
29  _trueOnly(iConfig.getUntrackedParameter<bool>("trueOnly", true)),
30  _runHisto(iConfig.getUntrackedParameter<bool>("runHisto", true)),
31  _runHistoProfile(iConfig.getUntrackedParameter<bool>("runHistoProfile", true)),
32  _runHistoBXProfile(iConfig.getUntrackedParameter<bool>("runHistoBXProfile", true)),
33  _runHistoBX2D(iConfig.getUntrackedParameter<bool>("runHistoBX2D", false)),
34  _histoParameters(iConfig.getUntrackedParameter<edm::ParameterSet>("histoParameters", edm::ParameterSet())),
35  _rhm(iC) {}
36 
38 
41  TFileDirectory* currdir = &(tfserv->tFileDirectory());
42 
43  if (!dirname.empty()) {
44  currdir = new TFileDirectory(tfserv->mkdir(dirname));
45  _currdir = currdir;
46  }
47 
48  edm::LogInfo("HistogramBooking") << "Vertex histogram booking in directory " << dirname;
49 
50  _hdeltax = currdir->make<TH1F>("deltax",
51  "(PV-BS) X position",
52  _histoParameters.getUntrackedParameter<unsigned int>("nBinX", 200),
53  _histoParameters.getUntrackedParameter<double>("xMin", -1.),
54  _histoParameters.getUntrackedParameter<double>("xMax", 1.));
55  _hdeltax->GetXaxis()->SetTitle("#Delta(X) [cm]");
56  _hdeltax->GetYaxis()->SetTitle("Vertices");
57 
58  _hdeltay = currdir->make<TH1F>("deltay",
59  "(PV-BS) Y position",
60  _histoParameters.getUntrackedParameter<unsigned int>("nBinY", 200),
61  _histoParameters.getUntrackedParameter<double>("yMin", -1.),
62  _histoParameters.getUntrackedParameter<double>("yMax", 1.));
63  _hdeltay->GetXaxis()->SetTitle("#Delta(Y) [cm]");
64  _hdeltay->GetYaxis()->SetTitle("Vertices");
65 
66  _hdeltaz = currdir->make<TH1F>("deltaz",
67  "(PV-BS) Z position",
68  _histoParameters.getUntrackedParameter<unsigned int>("nBinZ", 200),
69  _histoParameters.getUntrackedParameter<double>("zMin", -20.),
70  _histoParameters.getUntrackedParameter<double>("zMax", 20.));
71  _hdeltaz->GetXaxis()->SetTitle("#Delta(Z) [cm]");
72  _hdeltaz->GetYaxis()->SetTitle("Vertices");
73 
74  _hdeltaxvsz = currdir->make<TProfile>("deltaxvsz",
75  "(PV-BS) X position vs Z",
76  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile", 40),
77  _histoParameters.getUntrackedParameter<double>("zMinProfile", -20.),
78  _histoParameters.getUntrackedParameter<double>("zMaxProfile", 20.));
79  _hdeltaxvsz->GetXaxis()->SetTitle("Z [cm]");
80  _hdeltaxvsz->GetYaxis()->SetTitle("#Delta(X) [cm]");
81 
82  _hdeltayvsz = currdir->make<TProfile>("deltayvsz",
83  "(PV-BS) Y position vs Z",
84  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile", 40),
85  _histoParameters.getUntrackedParameter<double>("zMinProfile", -20.),
86  _histoParameters.getUntrackedParameter<double>("zMaxProfile", 20.));
87  _hdeltayvsz->GetXaxis()->SetTitle("Z [cm]");
88  _hdeltayvsz->GetYaxis()->SetTitle("#Delta(Y) [cm]");
89 
90  if (_runHisto) {
91  _hdeltaxrun = _rhm.makeTH1F("deltaxrun",
92  "(PV-BS) X position",
93  _histoParameters.getUntrackedParameter<unsigned int>("nBinX", 200),
94  _histoParameters.getUntrackedParameter<double>("xMin", -1.),
95  _histoParameters.getUntrackedParameter<double>("xMax", 1.));
96 
97  _hdeltayrun = _rhm.makeTH1F("deltayrun",
98  "(PV-BS) Y position",
99  _histoParameters.getUntrackedParameter<unsigned int>("nBinY", 200),
100  _histoParameters.getUntrackedParameter<double>("yMin", -1.),
101  _histoParameters.getUntrackedParameter<double>("yMax", 1.));
102 
103  _hdeltazrun = _rhm.makeTH1F("deltazrun",
104  "(PV-BS) Z position",
105  _histoParameters.getUntrackedParameter<unsigned int>("nBinZ", 200),
106  _histoParameters.getUntrackedParameter<double>("zMin", -20.),
107  _histoParameters.getUntrackedParameter<double>("zMax", 20.));
108 
109  _hdeltaxvszrun = _rhm.makeTProfile("deltaxvszrun",
110  "(PV-BS) X position vs Z",
111  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile", 40),
112  _histoParameters.getUntrackedParameter<double>("zMinProfile", -20.),
113  _histoParameters.getUntrackedParameter<double>("zMaxProfile", 20.));
114 
115  _hdeltayvszrun = _rhm.makeTProfile("deltayvszrun",
116  "(PV-BS) Y position vs Z",
117  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile", 40),
118  _histoParameters.getUntrackedParameter<double>("zMinProfile", -20.),
119  _histoParameters.getUntrackedParameter<double>("zMaxProfile", 20.));
120 
121  if (_runHistoProfile) {
123  "deltaxvsorbrun", "(PV-BS) X position vs orbit number", 4 * m_maxLS, 0.5, m_maxLS * 262144 + 0.5);
125  "deltayvsorbrun", "(PV-BS) Y position vs orbit number", 4 * m_maxLS, 0.5, m_maxLS * 262144 + 0.5);
127  "deltazvsorbrun", "(PV-BS) Z position vs orbit number", 4 * m_maxLS, 0.5, m_maxLS * 262144 + 0.5);
128  }
129  if (_runHistoBXProfile) {
130  _hdeltaxvsbxrun = _rhm.makeTProfile("deltaxvsbxrun", "(PV-BS) X position vs BX number", 3564, -0.5, 3563.5);
131  _hdeltayvsbxrun = _rhm.makeTProfile("deltayvsbxrun", "(PV-BS) Y position vs BX number", 3564, -0.5, 3563.5);
132  _hdeltazvsbxrun = _rhm.makeTProfile("deltazvsbxrun", "(PV-BS) Z position vs BX number", 3564, -0.5, 3563.5);
133  if (_runHistoBX2D) {
134  _hdeltaxvsbx2drun = _rhm.makeTH2F("deltaxvsbx2drun",
135  "(PV-BS) X position vs BX number",
136  3564,
137  -0.5,
138  3563.5,
139  _histoParameters.getUntrackedParameter<unsigned int>("nBinX", 200),
140  _histoParameters.getUntrackedParameter<double>("xMin", -1.),
141  _histoParameters.getUntrackedParameter<double>("xMax", 1.));
142  _hdeltayvsbx2drun = _rhm.makeTH2F("deltayvsbx2drun",
143  "(PV-BS) Y position vs BX number",
144  3564,
145  -0.5,
146  3563.5,
147  _histoParameters.getUntrackedParameter<unsigned int>("nBinY", 200),
148  _histoParameters.getUntrackedParameter<double>("yMin", -1.),
149  _histoParameters.getUntrackedParameter<double>("yMax", 1.));
150  _hdeltazvsbx2drun = _rhm.makeTH2F("deltazvsbx2drun",
151  "(PV-BS) Z position vs BX number",
152  3564,
153  -0.5,
154  3563.5,
155  _histoParameters.getUntrackedParameter<unsigned int>("nBinZ", 200),
156  _histoParameters.getUntrackedParameter<double>("zMin", -20.),
157  _histoParameters.getUntrackedParameter<double>("zMax", 20.));
158  }
159  }
160  }
161 }
162 
163 void BSvsPVHistogramMaker::beginRun(const unsigned int nrun) {
164  char runname[100];
165  sprintf(runname, "run_%d", nrun);
166 
167  TFileDirectory* currdir = _currdir;
168  if (currdir == nullptr) {
170  currdir = &(tfserv->tFileDirectory());
171  }
172 
173  _rhm.beginRun(nrun, *currdir);
174 
175  if (_runHisto) {
176  (*_hdeltaxrun)->GetXaxis()->SetTitle("#Delta(X) [cm]");
177  (*_hdeltaxrun)->GetYaxis()->SetTitle("Vertices");
178  (*_hdeltayrun)->GetXaxis()->SetTitle("#Delta(Y) [cm]");
179  (*_hdeltayrun)->GetYaxis()->SetTitle("Vertices");
180  (*_hdeltazrun)->GetXaxis()->SetTitle("#Delta(Z) [cm]");
181  (*_hdeltazrun)->GetYaxis()->SetTitle("Vertices");
182  (*_hdeltaxvszrun)->GetXaxis()->SetTitle("Z [cm]");
183  (*_hdeltaxvszrun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
184  (*_hdeltayvszrun)->GetXaxis()->SetTitle("Z [cm]");
185  (*_hdeltayvszrun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
186 
187  if (_runHistoProfile) {
188  (*_hdeltaxvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
189  (*_hdeltaxvsorbrun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
190  (*_hdeltaxvsorbrun)->SetCanExtend(TH1::kAllAxes);
191  (*_hdeltayvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
192  (*_hdeltayvsorbrun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
193  (*_hdeltayvsorbrun)->SetCanExtend(TH1::kAllAxes);
194  (*_hdeltazvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
195  (*_hdeltazvsorbrun)->GetYaxis()->SetTitle("#Delta(Z) [cm]");
196  (*_hdeltazvsorbrun)->SetCanExtend(TH1::kAllAxes);
197  }
198  if (_runHistoBXProfile) {
199  (*_hdeltaxvsbxrun)->GetXaxis()->SetTitle("BX");
200  (*_hdeltaxvsbxrun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
201  (*_hdeltayvsbxrun)->GetXaxis()->SetTitle("BX");
202  (*_hdeltayvsbxrun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
203  (*_hdeltazvsbxrun)->GetXaxis()->SetTitle("BX");
204  (*_hdeltazvsbxrun)->GetYaxis()->SetTitle("#Delta(Z) [cm]");
205  if (_runHistoBX2D) {
206  (*_hdeltaxvsbx2drun)->GetXaxis()->SetTitle("BX");
207  (*_hdeltaxvsbx2drun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
208  (*_hdeltayvsbx2drun)->GetXaxis()->SetTitle("BX");
209  (*_hdeltayvsbx2drun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
210  (*_hdeltazvsbx2drun)->GetXaxis()->SetTitle("BX");
211  (*_hdeltazvsbx2drun)->GetYaxis()->SetTitle("#Delta(Z) [cm]");
212  }
213  }
214  }
215 }
216 
217 void BSvsPVHistogramMaker::fill(const unsigned int orbit,
218  const int bx,
220  const reco::BeamSpot& bs) {
221  for (reco::VertexCollection::const_iterator vtx = vertices.begin(); vtx != vertices.end(); ++vtx) {
222  if (!(_trueOnly && vtx->isFake())) {
223  /*
224  double deltax = vtx->x()-bs.x0();
225  double deltay = vtx->y()-bs.y0();
226  double deltaz = vtx->z()-bs.z0();
227  */
228  double deltax = vtx->x() - x(bs, vtx->z());
229  double deltay = vtx->y() - y(bs, vtx->z());
230  double deltaz = vtx->z() - bs.z0();
231 
232  _hdeltax->Fill(deltax);
233  _hdeltay->Fill(deltay);
234  _hdeltaz->Fill(deltaz);
235  _hdeltaxvsz->Fill(vtx->z(), deltax);
236  _hdeltayvsz->Fill(vtx->z(), deltay);
237 
238  if (_runHisto) {
239  if (_hdeltaxrun && *_hdeltaxrun)
240  (*_hdeltaxrun)->Fill(deltax);
241  if (_hdeltayrun && *_hdeltayrun)
242  (*_hdeltayrun)->Fill(deltay);
243  if (_hdeltazrun && *_hdeltazrun)
244  (*_hdeltazrun)->Fill(deltaz);
246  (*_hdeltaxvszrun)->Fill(vtx->z(), deltax);
248  (*_hdeltayvszrun)->Fill(vtx->z(), deltay);
249  if (_runHistoProfile) {
251  (*_hdeltaxvsorbrun)->Fill(orbit, deltax);
253  (*_hdeltayvsorbrun)->Fill(orbit, deltay);
255  (*_hdeltazvsorbrun)->Fill(orbit, deltaz);
256  }
257  if (_runHistoBXProfile) {
259  (*_hdeltaxvsbxrun)->Fill(bx % 3564, deltax);
261  (*_hdeltayvsbxrun)->Fill(bx % 3564, deltay);
263  (*_hdeltazvsbxrun)->Fill(bx % 3564, deltaz);
264  if (_runHistoBX2D) {
266  (*_hdeltaxvsbx2drun)->Fill(bx % 3564, deltax);
268  (*_hdeltayvsbx2drun)->Fill(bx % 3564, deltay);
270  (*_hdeltazvsbx2drun)->Fill(bx % 3564, deltaz);
271  }
272  }
273  }
274  }
275  }
276 }
277 
280  const reco::BeamSpot& bs) {
281  fill(iEvent.orbitNumber(), iEvent.bunchCrossing(), vertices, bs);
282 }
283 
284 double BSvsPVHistogramMaker::x(const reco::BeamSpot& bs, const double z) const {
285  double x = bs.x0();
286 
287  // if(useSlope_) x += bs.dxdz()*z;
288  if (useSlope_)
289  x += bs.dxdz() * (z - bs.z0());
290 
291  // if(useSlope_) x = bs.x(z);
292 
293  return x;
294 }
295 
296 double BSvsPVHistogramMaker::y(const reco::BeamSpot& bs, const double z) const {
297  double y = bs.y0();
298 
299  // if(useSlope_) y += bs.dydz()*z;
300  if (useSlope_)
301  y += bs.dydz() * (z - bs.z0());
302 
303  // if(useSlope_) y = bs.y(z);
304 
305  return y;
306 }
BSvsPVHistogramMaker(edm::ConsumesCollector &&iC)
double z0() const
z coordinate
Definition: BeamSpot.h:65
T getUntrackedParameter(std::string const &, T const &) const
TFileDirectory * _currdir
const unsigned int m_maxLS
const edm::ParameterSet _histoParameters
int bunchCrossing() const
Definition: EventBase.h:64
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
RunHistogramManager _rhm
float float float z
_rhm(consumesCollector())
double dydz() const
dydz slope
Definition: BeamSpot.h:80
TFileDirectory & tFileDirectory()
Definition: TFileService.h:42
TH1F ** makeTH1F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
int iEvent
Definition: GenABIO.cc:224
TProfile ** makeTProfile(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
T * make(const Args &...args) const
make new ROOT object
int orbitNumber() const
Definition: EventBase.h:65
double dxdz() const
dxdz slope
Definition: BeamSpot.h:78
void book(const std::string dirname="")
Log< level::Info, false > LogInfo
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
double y(const reco::BeamSpot &bs, const double z) const
void beginRun(const edm::Run &iRun)
double y0() const
y coordinate
Definition: BeamSpot.h:63
TH2F ** makeTH2F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax, const unsigned int nbiny, const double ymin, const double ymax)
void beginRun(const unsigned int nrun)
void fill(const unsigned int orbit, const int bx, const reco::VertexCollection &vertices, const reco::BeamSpot &bs)
double x(const reco::BeamSpot &bs, const double z) const
double x0() const
x coordinate
Definition: BeamSpot.h:61