CMS 3D CMS Logo

Functions
MultiplicityPlotMacros.cc File Reference
#include <iostream>
#include <algorithm>
#include <vector>
#include <string>
#include "MultiplicityPlotMacros.h"
#include "DPGAnalysis/SiStripTools/interface/CommonAnalyzer.h"
#include "TFile.h"
#include "TH1D.h"
#include "TProfile.h"
#include "TDirectory.h"
#include "TCanvas.h"
#include "TStyle.h"
#include "TLegend.h"

Go to the source code of this file.

Functions

TH1D * AverageRunMultiplicity (TFile &ff, const char *module, const bool excludeLastBins=false, const char *histo="nTKdigivsorbrun")
 
void PlotPixelMultVtxPos (TFile *ff, const char *module)
 

Function Documentation

TH1D* AverageRunMultiplicity ( TFile &  ff,
const char *  module,
const bool  excludeLastBins = false,
const char *  histo = "nTKdigivsorbrun" 
)

Definition at line 101 of file MultiplicityPlotMacros.cc.

References gather_cfg::cout, CommonAnalyzer::getObject(), CommonAnalyzer::getRunList(), trackerHits::histo, mps_fire::i, writedatasetfile::runs, and CommonAnalyzer::setPath().

101  {
102 
103  CommonAnalyzer camult(&ff,"",module);
104 
105  TH1D* clusmult = new TH1D("clusmult","Average Multiplicity vs run",10,0.,10.);
106  clusmult->SetCanExtend(TH1::kXaxis);
107 
108  std::vector<unsigned int> runs = camult.getRunList();
109  std::sort(runs.begin(),runs.end());
110 
111  {
112  for(unsigned int i=0;i<runs.size();++i) {
113 
114  char runlabel[100];
115  sprintf(runlabel,"%d",runs[i]);
116  char runpath[100];
117  sprintf(runpath,"run_%d",runs[i]);
118  camult.setPath(runpath);
119 
120 
121  TProfile* multvstime=nullptr;
122  if(multvstime==nullptr) multvstime = (TProfile*)camult.getObject(histo);
123  if(multvstime) {
124  // compute mean exlucing the last filled bins
125 
126  if(excludeLastBins) {
127  int lastbin= multvstime->GetNbinsX()+1;
128  int firstbin= 1;
129  for(int ibin=multvstime->GetNbinsX()+1;ibin>0;--ibin) {
130  if(multvstime->GetBinEntries(ibin)!=0) {
131  lastbin=ibin;
132  break;
133  }
134  }
135 
136  std::cout << "Restricted range: " << firstbin << " " << lastbin << std::endl;
137  multvstime->GetXaxis()->SetRangeUser(multvstime->GetBinLowEdge(firstbin),multvstime->GetBinLowEdge(lastbin-2));
138  }
139  // fill the summary
140  clusmult->Fill(runlabel,multvstime->GetMean(2));
141 
142  }
143  }
144  }
145  return clusmult;
146 }
Definition: vlib.h:208
void PlotPixelMultVtxPos ( TFile *  ff,
const char *  module 
)

Definition at line 15 of file MultiplicityPlotMacros.cc.

References CommonAnalyzer::getObject(), mps_fire::i, tablePrinter::labels, create_public_lumi_plots::leg, callgraph::path, CommonAnalyzer::setPath(), and AlCaHLTBitMon_QueryRunRegistry::string.

15  {
16 
17  CommonAnalyzer camult(ff,"",module);
18  // camult.setPath("VtxPosCorr");
19 
20  std::vector<std::string> labels;
21  labels.push_back("FPIX_m");
22  labels.push_back("BPIX_L1_mod_1");
23  labels.push_back("BPIX_L1_mod_2");
24  labels.push_back("BPIX_L1_mod_3");
25  labels.push_back("BPIX_L1_mod_4");
26  labels.push_back("BPIX_L1_mod_5");
27  labels.push_back("BPIX_L1_mod_6");
28  labels.push_back("BPIX_L1_mod_7");
29  labels.push_back("BPIX_L1_mod_8");
30  labels.push_back("FPIX_p");
31  labels.push_back("BPIX_L1");
32  labels.push_back("BPIX_L2");
33  labels.push_back("BPIX_L3");
34  labels.push_back("Lumi");
35 
36  std::vector<TProfile*> profs;
37 
38  for(unsigned int i=0;i<labels.size();++i) {
39 
40  std::string path = "VtxPosCorr/"+labels[i];
41  camult.setPath(path.c_str());
42 
43  std::string hname = "n"+labels[i]+"digivsvtxposprof";
44  profs.push_back((TProfile*)camult.getObject(hname.c_str()));
45 
46  }
47 
48  TCanvas* cc = new TCanvas("BPIX L1 details","BPIX L1 details",1000,1000);
49  gPad->Divide(2,2);
50 
51  for(unsigned int i = 1;i<5;++i) {
52  cc->cd(i);
53  if(profs[i] && profs[9-i]) {
54  profs[i]->Draw();
55  profs[9-i]->SetLineColor(kRed);
56  profs[9-i]->SetMarkerColor(kRed);
57  profs[9-i]->Draw("same");
58  TLegend* leg = new TLegend(0.4,0.8,0.6,0.9,"Occupancy");
59  leg->SetFillStyle(0);
60  leg->AddEntry(profs[i],labels[i].c_str(),"l");
61  leg->AddEntry(profs[9-i],labels[9-i].c_str(),"l");
62  leg->Draw();
63  }
64  }
65  new TCanvas("FPIX","FPIX");
66  if(profs[0] && profs[9]) {
67  profs[0]->Draw();
68  profs[9]->SetLineColor(kRed);
69  profs[9]->SetMarkerColor(kRed);
70  profs[9]->Draw("same");
71  TLegend* leg = new TLegend(0.4,0.8,0.6,0.9,"Occupancy");
72  leg->SetFillStyle(0);
73  leg->AddEntry(profs[0],labels[0].c_str(),"l");
74  leg->AddEntry(profs[9],labels[9].c_str(),"l");
75  leg->Draw();
76  }
77 
78  gStyle->SetOptStat(11);
79  gStyle->SetOptFit(11);
80  new TCanvas("BPIXL1","BPIX L1");
81  profs[10]->Fit("pol2");
82  new TCanvas("BPIXL2","BPIX L2");
83  profs[11]->Fit("pol2");
84  new TCanvas("BPIXL3","BPIX L3");
85  profs[12]->Fit("pol2");
86 
87  new TCanvas("LumiAdd","LumiAdd");
88  TH1D* hlumi = profs[11]->ProjectionX("lumi");
89  TH1D* hbpixl3 = profs[12]->ProjectionX("bpixl3");
90  TH1D* hfpixm = profs[0]->ProjectionX("fpixm");
91  TH1D* hfpixp = profs[9]->ProjectionX("fpixp");
92  hlumi->SetTitle("BPIX L2+L3 + FPIX multiplicity vs vtx z position");
93  hlumi->Add(hbpixl3);
94  hlumi->Add(hfpixm);
95  hlumi->Add(hfpixp);
96  hlumi->Fit("pol2");
97  new TCanvas("Lumi","Lumi");
98  profs[13]->Fit("pol2");
99 }
Definition: vlib.h:208