CMS 3D CMS Logo

BSvsPVHistogramMaker.cc
Go to the documentation of this file.
9 #include "TH1F.h"
10 #include "TH2F.h"
11 #include "TProfile.h"
12 
13 
15  _currdir(0), m_maxLS(100), useSlope_(true), _trueOnly(true),
16  _runHisto(true), _runHistoProfile(true), _runHistoBXProfile(true), _runHistoBX2D(false), _histoParameters(), _rhm(iC) { }
17 
19  _currdir(0),
20  m_maxLS(iConfig.getParameter<unsigned int>("maxLSBeforeRebin")),
21  useSlope_(iConfig.getParameter<bool>("useSlope")),
22  _trueOnly(iConfig.getUntrackedParameter<bool>("trueOnly",true)),
23  _runHisto(iConfig.getUntrackedParameter<bool>("runHisto",true)),
24  _runHistoProfile(iConfig.getUntrackedParameter<bool>("runHistoProfile",true)),
25  _runHistoBXProfile(iConfig.getUntrackedParameter<bool>("runHistoBXProfile",true)),
26  _runHistoBX2D(iConfig.getUntrackedParameter<bool>("runHistoBX2D",false)),
27  _histoParameters(iConfig.getUntrackedParameter<edm::ParameterSet>("histoParameters",edm::ParameterSet())),
28  _rhm(iC)
29 { }
30 
31 
33 
34  delete _currdir;
35 
36 }
37 
38 
40 
42  TFileDirectory* currdir = &(tfserv->tFileDirectory());
43 
44  if(dirname!="") {
45  currdir = new TFileDirectory(tfserv->mkdir(dirname));
46  _currdir = currdir;
47  }
48 
49  edm::LogInfo("HistogramBooking") << "Vertex histogram booking in directory " << dirname;
50 
51  _hdeltax = currdir->make<TH1F>("deltax","(PV-BS) X position",
52  _histoParameters.getUntrackedParameter<unsigned int>("nBinX",200),
53  _histoParameters.getUntrackedParameter<double>("xMin",-1.),
54  _histoParameters.getUntrackedParameter<double>("xMax",1.)
55  );
56  _hdeltax->GetXaxis()->SetTitle("#Delta(X) [cm]"); _hdeltax->GetYaxis()->SetTitle("Vertices");
57 
58  _hdeltay = currdir->make<TH1F>("deltay","(PV-BS) Y position",
59  _histoParameters.getUntrackedParameter<unsigned int>("nBinY",200),
60  _histoParameters.getUntrackedParameter<double>("yMin",-1.),
61  _histoParameters.getUntrackedParameter<double>("yMax",1.)
62  );
63  _hdeltay->GetXaxis()->SetTitle("#Delta(Y) [cm]"); _hdeltay->GetYaxis()->SetTitle("Vertices");
64 
65  _hdeltaz = currdir->make<TH1F>("deltaz","(PV-BS) Z position",
66  _histoParameters.getUntrackedParameter<unsigned int>("nBinZ",200),
67  _histoParameters.getUntrackedParameter<double>("zMin",-20.),
68  _histoParameters.getUntrackedParameter<double>("zMax",20.)
69  );
70  _hdeltaz->GetXaxis()->SetTitle("#Delta(Z) [cm]"); _hdeltaz->GetYaxis()->SetTitle("Vertices");
71 
72  _hdeltaxvsz = currdir->make<TProfile>("deltaxvsz","(PV-BS) X position vs Z",
73  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile",40),
74  _histoParameters.getUntrackedParameter<double>("zMinProfile",-20.),
75  _histoParameters.getUntrackedParameter<double>("zMaxProfile",20.)
76  );
77  _hdeltaxvsz->GetXaxis()->SetTitle("Z [cm]"); _hdeltaxvsz->GetYaxis()->SetTitle("#Delta(X) [cm]");
78 
79  _hdeltayvsz = currdir->make<TProfile>("deltayvsz","(PV-BS) Y position vs Z",
80  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile",40),
81  _histoParameters.getUntrackedParameter<double>("zMinProfile",-20.),
82  _histoParameters.getUntrackedParameter<double>("zMaxProfile",20.)
83  );
84  _hdeltayvsz->GetXaxis()->SetTitle("Z [cm]"); _hdeltayvsz->GetYaxis()->SetTitle("#Delta(Y) [cm]");
85 
86 
87 
88 
89  if(_runHisto) {
90  _hdeltaxrun = _rhm.makeTH1F("deltaxrun","(PV-BS) X position",
91  _histoParameters.getUntrackedParameter<unsigned int>("nBinX",200),
92  _histoParameters.getUntrackedParameter<double>("xMin",-1.),
93  _histoParameters.getUntrackedParameter<double>("xMax",1.));
94 
95  _hdeltayrun = _rhm.makeTH1F("deltayrun","(PV-BS) Y position",
96  _histoParameters.getUntrackedParameter<unsigned int>("nBinY",200),
97  _histoParameters.getUntrackedParameter<double>("yMin",-1.),
98  _histoParameters.getUntrackedParameter<double>("yMax",1.));
99 
100  _hdeltazrun = _rhm.makeTH1F("deltazrun","(PV-BS) Z position",
101  _histoParameters.getUntrackedParameter<unsigned int>("nBinZ",200),
102  _histoParameters.getUntrackedParameter<double>("zMin",-20.),
103  _histoParameters.getUntrackedParameter<double>("zMax",20.));
104 
105  _hdeltaxvszrun = _rhm.makeTProfile("deltaxvszrun","(PV-BS) X position vs Z",
106  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile",40),
107  _histoParameters.getUntrackedParameter<double>("zMinProfile",-20.),
108  _histoParameters.getUntrackedParameter<double>("zMaxProfile",20.)
109  );
110 
111  _hdeltayvszrun = _rhm.makeTProfile("deltayvszrun","(PV-BS) Y position vs Z",
112  _histoParameters.getUntrackedParameter<unsigned int>("nBinZProfile",40),
113  _histoParameters.getUntrackedParameter<double>("zMinProfile",-20.),
114  _histoParameters.getUntrackedParameter<double>("zMaxProfile",20.)
115  );
116 
117  if(_runHistoProfile) {
118  _hdeltaxvsorbrun = _rhm.makeTProfile("deltaxvsorbrun","(PV-BS) X position vs orbit number",4*m_maxLS,0.5,m_maxLS*262144+0.5);
119  _hdeltayvsorbrun = _rhm.makeTProfile("deltayvsorbrun","(PV-BS) Y position vs orbit number",4*m_maxLS,0.5,m_maxLS*262144+0.5);
120  _hdeltazvsorbrun = _rhm.makeTProfile("deltazvsorbrun","(PV-BS) Z position vs orbit number",4*m_maxLS,0.5,m_maxLS*262144+0.5);
121  }
122  if(_runHistoBXProfile) {
123  _hdeltaxvsbxrun = _rhm.makeTProfile("deltaxvsbxrun","(PV-BS) X position vs BX number",3564,-0.5,3563.5);
124  _hdeltayvsbxrun = _rhm.makeTProfile("deltayvsbxrun","(PV-BS) Y position vs BX number",3564,-0.5,3563.5);
125  _hdeltazvsbxrun = _rhm.makeTProfile("deltazvsbxrun","(PV-BS) Z position vs BX number",3564,-0.5,3563.5);
126  if(_runHistoBX2D) {
127  _hdeltaxvsbx2drun = _rhm.makeTH2F("deltaxvsbx2drun","(PV-BS) X position vs BX number",3564,-0.5,3563.5,
128  _histoParameters.getUntrackedParameter<unsigned int>("nBinX",200),
129  _histoParameters.getUntrackedParameter<double>("xMin",-1.),
130  _histoParameters.getUntrackedParameter<double>("xMax",1.));
131  _hdeltayvsbx2drun = _rhm.makeTH2F("deltayvsbx2drun","(PV-BS) Y position vs BX number",3564,-0.5,3563.5,
132  _histoParameters.getUntrackedParameter<unsigned int>("nBinY",200),
133  _histoParameters.getUntrackedParameter<double>("yMin",-1.),
134  _histoParameters.getUntrackedParameter<double>("yMax",1.));
135  _hdeltazvsbx2drun = _rhm.makeTH2F("deltazvsbx2drun","(PV-BS) Z position vs BX number",3564,-0.5,3563.5,
136  _histoParameters.getUntrackedParameter<unsigned int>("nBinZ",200),
137  _histoParameters.getUntrackedParameter<double>("zMin",-20.),
138  _histoParameters.getUntrackedParameter<double>("zMax",20.));
139  }
140  }
141 
142  }
143 }
144 
145 void BSvsPVHistogramMaker::beginRun(const unsigned int nrun) {
146 
147  char runname[100];
148  sprintf(runname,"run_%d",nrun);
149 
150  TFileDirectory* currdir = _currdir;
151  if(currdir==0) {
153  currdir = &(tfserv->tFileDirectory());
154  }
155 
156  _rhm.beginRun(nrun,*currdir);
157 
158  if(_runHisto) {
159  (*_hdeltaxrun)->GetXaxis()->SetTitle("#Delta(X) [cm]"); (*_hdeltaxrun)->GetYaxis()->SetTitle("Vertices");
160  (*_hdeltayrun)->GetXaxis()->SetTitle("#Delta(Y) [cm]"); (*_hdeltayrun)->GetYaxis()->SetTitle("Vertices");
161  (*_hdeltazrun)->GetXaxis()->SetTitle("#Delta(Z) [cm]"); (*_hdeltazrun)->GetYaxis()->SetTitle("Vertices");
162  (*_hdeltaxvszrun)->GetXaxis()->SetTitle("Z [cm]"); (*_hdeltaxvszrun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
163  (*_hdeltayvszrun)->GetXaxis()->SetTitle("Z [cm]"); (*_hdeltayvszrun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
164 
165  if(_runHistoProfile) {
166  (*_hdeltaxvsorbrun)->GetXaxis()->SetTitle("time [orbit#]"); (*_hdeltaxvsorbrun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
167  (*_hdeltaxvsorbrun)->SetCanExtend(TH1::kAllAxes);
168  (*_hdeltayvsorbrun)->GetXaxis()->SetTitle("time [orbit#]"); (*_hdeltayvsorbrun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
169  (*_hdeltayvsorbrun)->SetCanExtend(TH1::kAllAxes);
170  (*_hdeltazvsorbrun)->GetXaxis()->SetTitle("time [orbit#]"); (*_hdeltazvsorbrun)->GetYaxis()->SetTitle("#Delta(Z) [cm]");
171  (*_hdeltazvsorbrun)->SetCanExtend(TH1::kAllAxes);
172  }
173  if(_runHistoBXProfile) {
174  (*_hdeltaxvsbxrun)->GetXaxis()->SetTitle("BX"); (*_hdeltaxvsbxrun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
175  (*_hdeltayvsbxrun)->GetXaxis()->SetTitle("BX"); (*_hdeltayvsbxrun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
176  (*_hdeltazvsbxrun)->GetXaxis()->SetTitle("BX"); (*_hdeltazvsbxrun)->GetYaxis()->SetTitle("#Delta(Z) [cm]");
177  if(_runHistoBX2D) {
178  (*_hdeltaxvsbx2drun)->GetXaxis()->SetTitle("BX"); (*_hdeltaxvsbx2drun)->GetYaxis()->SetTitle("#Delta(X) [cm]");
179  (*_hdeltayvsbx2drun)->GetXaxis()->SetTitle("BX"); (*_hdeltayvsbx2drun)->GetYaxis()->SetTitle("#Delta(Y) [cm]");
180  (*_hdeltazvsbx2drun)->GetXaxis()->SetTitle("BX"); (*_hdeltazvsbx2drun)->GetYaxis()->SetTitle("#Delta(Z) [cm]");
181  }
182  }
183 
184  }
185 }
186 
187 void BSvsPVHistogramMaker::fill(const unsigned int orbit, const int bx, const reco::VertexCollection& vertices, const reco::BeamSpot& bs) {
188 
189  for(reco::VertexCollection::const_iterator vtx=vertices.begin();vtx!=vertices.end();++vtx) {
190 
191  if(!(_trueOnly && vtx->isFake())) {
192 
193  /*
194  double deltax = vtx->x()-bs.x0();
195  double deltay = vtx->y()-bs.y0();
196  double deltaz = vtx->z()-bs.z0();
197  */
198  double deltax = vtx->x()-x(bs,vtx->z());
199  double deltay = vtx->y()-y(bs,vtx->z());
200  double deltaz = vtx->z()-bs.z0();
201 
202  _hdeltax->Fill(deltax);
203  _hdeltay->Fill(deltay);
204  _hdeltaz->Fill(deltaz);
205  _hdeltaxvsz->Fill(vtx->z(),deltax);
206  _hdeltayvsz->Fill(vtx->z(),deltay);
207 
208  if(_runHisto) {
209  if(_hdeltaxrun && *_hdeltaxrun ) (*_hdeltaxrun)->Fill(deltax);
210  if(_hdeltayrun && *_hdeltayrun ) (*_hdeltayrun)->Fill(deltay);
211  if(_hdeltazrun && *_hdeltazrun ) (*_hdeltazrun)->Fill(deltaz);
212  if(_hdeltaxvszrun && *_hdeltaxvszrun ) (*_hdeltaxvszrun)->Fill(vtx->z(),deltax);
213  if(_hdeltayvszrun && *_hdeltayvszrun ) (*_hdeltayvszrun)->Fill(vtx->z(),deltay);
214  if(_runHistoProfile) {
215  if(_hdeltaxvsorbrun && *_hdeltaxvsorbrun ) (*_hdeltaxvsorbrun)->Fill(orbit,deltax);
216  if(_hdeltayvsorbrun && *_hdeltayvsorbrun ) (*_hdeltayvsorbrun)->Fill(orbit,deltay);
217  if(_hdeltazvsorbrun && *_hdeltazvsorbrun ) (*_hdeltazvsorbrun)->Fill(orbit,deltaz);
218  }
219  if(_runHistoBXProfile) {
220  if(_hdeltaxvsbxrun && *_hdeltaxvsbxrun ) (*_hdeltaxvsbxrun)->Fill(bx%3564,deltax);
221  if(_hdeltayvsbxrun && *_hdeltayvsbxrun ) (*_hdeltayvsbxrun)->Fill(bx%3564,deltay);
222  if(_hdeltazvsbxrun && *_hdeltazvsbxrun ) (*_hdeltazvsbxrun)->Fill(bx%3564,deltaz);
223  if(_runHistoBX2D) {
224  if(_hdeltaxvsbx2drun && *_hdeltaxvsbx2drun ) (*_hdeltaxvsbx2drun)->Fill(bx%3564,deltax);
225  if(_hdeltayvsbx2drun && *_hdeltayvsbx2drun ) (*_hdeltayvsbx2drun)->Fill(bx%3564,deltay);
226  if(_hdeltazvsbx2drun && *_hdeltazvsbx2drun ) (*_hdeltazvsbx2drun)->Fill(bx%3564,deltaz);
227  }
228  }
229  }
230  }
231  }
232 }
233 
235 
236  fill(iEvent.orbitNumber(),iEvent.bunchCrossing(),vertices,bs);
237 
238 }
239 
240 double BSvsPVHistogramMaker::x(const reco::BeamSpot& bs, const double z) const {
241 
242  double x = bs.x0();
243 
244  // if(useSlope_) x += bs.dxdz()*z;
245  if(useSlope_) x += bs.dxdz()*(z-bs.z0());
246 
247  // if(useSlope_) x = bs.x(z);
248 
249  return x;
250 
251 }
252 
253 double BSvsPVHistogramMaker::y(const reco::BeamSpot& bs, const double z) const {
254 
255  double y = bs.y0();
256 
257  // if(useSlope_) y += bs.dydz()*z;
258  if(useSlope_) y += bs.dydz()*(z-bs.z0());
259 
260  // if(useSlope_) y = bs.y(z);
261 
262  return y;
263 
264 }
265 
BSvsPVHistogramMaker(edm::ConsumesCollector &&iC)
double z0() const
z coordinate
Definition: BeamSpot.h:68
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:66
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
RunHistogramManager _rhm
double dydz() const
dydz slope
Definition: BeamSpot.h:84
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:230
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:67
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
void book(const std::string dirname="")
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)
HLT enums.
double y0() const
y coordinate
Definition: BeamSpot.h:66
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:64