CMS 3D CMS Logo

BSvsBPIX.cc
Go to the documentation of this file.
2 #include "TFile.h"
3 #include "TH1F.h"
4 #include "TProfile.h"
5 #include "TGraphErrors.h"
6 #include "TCanvas.h"
7 #include "TDirectory.h"
8 #include <iostream>
9 
10 void BSvsBPIXPlot(TFile* ff, const char* bsmodule, const char* occumodule, const int run) {
11  TGraphErrors* bspos = new TGraphErrors();
12  TGraphErrors* bpixpos = new TGraphErrors();
13 
14  if (ff) {
15  char bsfolder[200];
16  sprintf(bsfolder, "%s/run_%d", bsmodule, run);
17  if (ff->cd(bsfolder)) {
18  TH1F* bsx = (TH1F*)gDirectory->Get("bsxrun");
19  TH1F* bsy = (TH1F*)gDirectory->Get("bsyrun");
20  if (bsx && bsy) {
21  std::cout << "beam spot position (" << bsx->GetMean() << "+/-" << bsx->GetMeanError() << "," << bsy->GetMean()
22  << "+/-" << bsy->GetMeanError() << ")" << std::endl;
23  bspos->SetPoint(0, bsx->GetMean(), bsy->GetMean());
24  bspos->SetPointError(0, bsx->GetMeanError(), bsy->GetMeanError());
25  }
26  }
27  char occufolder[200];
28  sprintf(occufolder, "%s/run_%d", occumodule, run);
29  if (ff->cd(occufolder)) {
30  TProfile* xmean = (TProfile*)gDirectory->Get("avex");
31  TProfile* ymean = (TProfile*)gDirectory->Get("avey");
32  if (xmean && ymean) {
33  for (int i = 1; i <= xmean->GetNbinsX(); ++i) {
34  if (xmean->GetBinEntries(i) > 0) {
35  std::cout << "ladder position " << i << " : (" << xmean->GetBinContent(i) << "+/-" << xmean->GetBinError(i)
36  << "," << ymean->GetBinContent(i) << "+/-" << ymean->GetBinError(i) << ")" << std::endl;
37  int point = bpixpos->GetN();
38  bpixpos->SetPoint(point, xmean->GetBinContent(i), ymean->GetBinContent(i));
39  bpixpos->SetPointError(point, xmean->GetBinError(i), ymean->GetBinError(i));
40  }
41  }
42  }
43  }
44  }
45  new TCanvas("bsbpix", "bsbpix", 500, 500);
46  bpixpos->Draw("ap");
47  bspos->Draw("p");
48 }
mps_fire.i
i
Definition: mps_fire.py:355
gather_cfg.cout
cout
Definition: gather_cfg.py:144
BSvsBPIXPlot
void BSvsBPIXPlot(TFile *ff, const char *bsmodule, const char *occumodule, const int run)
Definition: BSvsBPIX.cc:10
alignCSCRings.ff
ff
Definition: alignCSCRings.py:148
writedatasetfile.run
run
Definition: writedatasetfile.py:27
BSvsBPIX.h
point
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5