CMS 3D CMS Logo

Functions
multibsvspvplots.cc File Reference
#include "multibsvspvplots.h"
#include <vector>
#include <algorithm>
#include <string>
#include <map>
#include <iostream>
#include "TPad.h"
#include "TH1D.h"
#include "TFile.h"
#include "TH2F.h"
#include "TH1F.h"
#include "TF1.h"
#include "TProfile.h"
#include "TProfile2D.h"
#include "DPGAnalysis/SiStripTools/interface/CommonAnalyzer.h"
#include "TCanvas.h"
#include "TSystem.h"
#include "TStyle.h"
#include "TText.h"
#include "TLegend.h"

Go to the source code of this file.

Functions

void multibsvspvplots (const char *module, const char *label, const char *postfix, const char *)
 

Function Documentation

void multibsvspvplots ( const char *  module,
const char *  label,
const char *  postfix,
const char *   
)

Definition at line 22 of file multibsvspvplots.cc.

References stringResolutionProvider_cfi::bin, gather_cfg::cout, FrontierConditions_GlobalTag_cff::file, CommonAnalyzer::getObject(), CommonAnalyzer::getRunList(), mps_fire::i, compare_using_db::ifile, gather_cfg::runs, CommonAnalyzer::setPath(), and AlCaHLTBitMon_QueryRunRegistry::string.

22  {
23 
24  TFile* file[5];
25 
26  file[0] = new TFile("rootfiles/Tracking_PFG_Run2012A_prompt_minbias_v1_190456-194076_muon_relumi_v55_fittedV0.root");
27  file[1] = new TFile("rootfiles/Tracking_PFG_Run2012B_prompt_minbias_v1_190456-196531_muon_relumi_v55_fittedV0.root");
28  file[2] = new TFile("rootfiles/Tracking_PFG_Run2012C_prompt_minbias_v1_190456-199011_muon_relumi_v55_fittedV0.root");
29  file[3] = new TFile("rootfiles/Tracking_PFG_Run2012C_prompt_minbias_v2_190456-203002_muon_relumi_v57_fittedV0.root");
30  file[4] = new TFile("rootfiles/Tracking_PFG_Run2012D_prompt_minbias_v1_190456-208686_muon_relumi_v57_fittedV0.root");
31 
32  char modfull[300];
33  sprintf(modfull,"%s%s",module,postfix);
34 
35  char labfull[300];
36  sprintf(labfull,"%s%s",label,postfix);
37 
38  // Summary histograms
39  TH1D* deltaxsum = new TH1D("deltaxsum","(PV-BS) Fitted X position vs run",10,0.,10.);
40  deltaxsum->SetCanExtend(TH1::kAllAxes);
41  TH1D* deltaysum = new TH1D("deltaysum","(PV-BS) Fitted Y position vs run",10,0.,10.);
42  deltaysum->SetCanExtend(TH1::kAllAxes);
43 
44  TH1D* deltaxmeansum = new TH1D("deltaxmeansum","(PV-BS) Mean X position vs run",10,0.,10.);
45  deltaxmeansum->SetCanExtend(TH1::kAllAxes);
46  TH1D* deltaymeansum = new TH1D("deltaymeansum","(PV-BS) Mean Y position vs run",10,0.,10.);
47  deltaymeansum->SetCanExtend(TH1::kAllAxes);
48 
49  TH1D* deltazmeansum = new TH1D("deltazmeansum","(PV-BS) Mean Z position vs run",10,0.,10.);
50  deltazmeansum->SetCanExtend(TH1::kAllAxes);
51 
52  TF1* fdoubleg = new TF1("doubleg","[1]*exp(-0.5*((x-[0])/[2])**2)+[3]*exp(-0.5*((x-[0])/[4])**2)+[5]*exp(sqrt((x-[0])**2)*[6])",-.2,.2);
53  fdoubleg->SetLineColor(kRed);
54  fdoubleg->SetLineWidth(1);
55 
56  for(unsigned int ifile=0;ifile<5;++ifile) {
57 
58  CommonAnalyzer castat(file[ifile],"",modfull);
59 
60  std::vector<unsigned int> runs = castat.getRunList();
61  std::sort(runs.begin(),runs.end());
62 
63  {
64 
65  std::cout << "Found " << runs.size() << " runs" << std::endl;
66 
67  for(unsigned int i=0;i<runs.size();++i) {
68 
69  char runlabel[100];
70  sprintf(runlabel,"%d",runs[i]);
71  char runpath[100];
72  sprintf(runpath,"run_%d",runs[i]);
73  castat.setPath(runpath);
74  std::cout << runpath << std::endl;
75 
76  TH1F* deltax = (TH1F*)castat.getObject("deltaxrun");
77  if (deltax && deltax->GetEntries()>0) {
78 
79  fdoubleg->SetParameter(0,deltax->GetMean());
80  fdoubleg->SetParameter(2,deltax->GetRMS());
81  fdoubleg->SetParameter(4,deltax->GetRMS());
82  fdoubleg->SetParameter(1,deltax->GetMaximum());
83  fdoubleg->SetParameter(3,0.1*deltax->GetMaximum());
84  fdoubleg->SetParameter(5,0.1*deltax->GetMaximum());
85  /* const int result = */ deltax->Fit(fdoubleg,"q0b","",-.05,.05);
86 
87  int bin = deltaxsum->Fill(runlabel,fdoubleg->GetParameter(0));
88  deltaxsum->SetBinError(bin,fdoubleg->GetParError(0));
89 
90  bin = deltaxmeansum->Fill(runlabel,deltax->GetMean());
91  deltaxmeansum->SetBinError(bin,deltax->GetMeanError());
92 
93  }
94  TH1F* deltay = (TH1F*)castat.getObject("deltayrun");
95  if (deltay && deltay->GetEntries()>0) {
96 
97  fdoubleg->SetParameter(0,deltay->GetMean());
98  fdoubleg->SetParameter(2,deltay->GetRMS());
99  fdoubleg->SetParameter(4,deltay->GetRMS());
100  fdoubleg->SetParameter(1,deltay->GetMaximum());
101  fdoubleg->SetParameter(3,0.1*deltay->GetMaximum());
102  fdoubleg->SetParameter(5,0.1*deltay->GetMaximum());
103  /* const int result = */ deltay->Fit(fdoubleg,"q0b","",-.05,.05);
104  int bin = deltaysum->Fill(runlabel,fdoubleg->GetParameter(0));
105  deltaysum->SetBinError(bin,fdoubleg->GetParError(0));
106 
107  bin = deltaymeansum->Fill(runlabel,deltay->GetMean());
108  deltaymeansum->SetBinError(bin,deltay->GetMeanError());
109 
110  }
111  TH1F* deltaz = (TH1F*)castat.getObject("deltazrun");
112  if (deltaz && deltaz->GetEntries()>0) {
113 
114  int bin = deltazmeansum->Fill(runlabel,deltaz->GetMean());
115  deltazmeansum->SetBinError(bin,deltaz->GetMeanError());
116 
117  }
118 
119 
120 
121  }
122 
123  }
124  }
125 
126  std::string plotfilename;
127 
128  plotfilename = "/afs/cern.ch/cms/tracking/output/";
129  // plotfilename += dirname;
130  // plotfilename += "/";
131  // plotfilename += labfull;
132  plotfilename += "/deltaxsum_";
133  plotfilename += labfull;
134  // plotfilename += "_";
135  // plotfilename += dirname;
136  plotfilename += ".gif";
137 
138  /* TCanvas * cwidedeltax = */ new TCanvas(plotfilename.c_str(),plotfilename.c_str(),1500,500);
139 
140  deltaxsum->SetLineColor(kRed); deltaxsum->SetMarkerColor(kRed);
141  deltaxsum->GetYaxis()->SetRangeUser(-.002,.002);
142  deltaxsum->GetYaxis()->SetTitle("#Delta x (cm)");
143  deltaxsum->Draw();
144  deltaxmeansum->Draw("esame");
145  TLegend deltaxleg(.7,.8,.85,.9,"#Delta(x)");
146  deltaxleg.AddEntry(deltaxsum,"fitted mean","l");
147  deltaxleg.AddEntry(deltaxmeansum,"aritm. mean","l");
148  deltaxleg.Draw();
149  gPad->Print(plotfilename.c_str());
150  // delete cwidedeltax;
151 
152  plotfilename = "/afs/cern.ch/cms/tracking/output/";
153  // plotfilename += dirname;
154  // plotfilename += "/";
155  // plotfilename += labfull;
156  plotfilename += "/deltaysum_";
157  plotfilename += labfull;
158  // plotfilename += "_";
159  // plotfilename += dirname;
160  plotfilename += ".gif";
161 
162  /* TCanvas * cwidedeltay = */ new TCanvas(plotfilename.c_str(),plotfilename.c_str(),1500,500);
163 
164  deltaysum->SetLineColor(kRed); deltaysum->SetMarkerColor(kRed);
165  deltaysum->GetYaxis()->SetRangeUser(-.002,.002);
166  deltaysum->GetYaxis()->SetTitle("#Delta y (cm)");
167  deltaysum->Draw();
168  deltaymeansum->Draw("esame");
169  TLegend deltayleg(.7,.8,.85,.9,"#Delta(y)");
170  deltayleg.AddEntry(deltaysum,"fitted mean","l");
171  deltayleg.AddEntry(deltaymeansum,"aritm. mean","l");
172  deltayleg.Draw();
173  gPad->Print(plotfilename.c_str());
174  // delete cwidedeltay;
175 
176 
177  plotfilename = "/afs/cern.ch/cms/tracking/output/";
178  // plotfilename += dirname;
179  // plotfilename += "/";
180  // plotfilename += labfull;
181  plotfilename += "/deltazsum_";
182  plotfilename += labfull;
183  // plotfilename += "_";
184  // plotfilename += dirname;
185  plotfilename += ".gif";
186 
187  /* TCanvas * cwidedeltaz = */ new TCanvas(plotfilename.c_str(),plotfilename.c_str(),1500,500);
188 
189  deltazmeansum->GetYaxis()->SetRangeUser(-2.,2.);
190  deltazmeansum->GetYaxis()->SetTitle("#Delta z (cm)");
191  deltazmeansum->Draw();
192  gPad->Print(plotfilename.c_str());
193  // delete cwidedeltaz;
194  /*
195  delete deltaxsum;
196  delete deltaysum;
197  delete deltaxmeansum;
198  delete deltaymeansum;
199  delete deltazmeansum;
200  */
201  delete fdoubleg;
202 
203 }
char const * label
bin
set the eta bin as selection string.
Definition: vlib.h:208