CMS 3D CMS Logo

Classes | Functions
OOTMultiplicityPlotMacros.h File Reference
#include "TH1F.h"

Go to the source code of this file.

Classes

class  OOTResult
 
class  OOTSummary
 

Functions

OOTResultComputeOOTFraction (TFile *ff, const char *itmodule, const char *ootmodule, const char *etmodule, const int run, const char *hname, const bool &perFill=false)
 
OOTSummaryComputeOOTFractionvsFill (TFile *ff, const char *itmodule, const char *ootmodule, const char *etmodule, const char *hname, OOTSummary *ootsumm=nullptr)
 
OOTSummaryComputeOOTFractionvsRun (TFile *ff, const char *itmodule, const char *ootmodule, const char *etmodule, const char *hname, OOTSummary *ootsumm=nullptr)
 
std::vector< int > FillingScheme (TFile *ff, const char *path, const float thr=0.)
 
std::vector< int > FillingSchemeFromProfile (TFile *ff, const char *path, const char *hname, const float thr=0.)
 

Function Documentation

OOTResult* ComputeOOTFraction ( TFile *  ff,
const char *  itmodule,
const char *  ootmodule,
const char *  etmodule,
const int  run,
const char *  hname,
const bool &  perFill = false 
)

Definition at line 75 of file OOTMultiplicityPlotMacros.cc.

References gather_cfg::cout, FillingScheme(), FillingSchemeFromProfile(), OOTResult::hratio, OOTResult::nfilledbx, OOTResult::ngoodbx, OOTResult::ootfrac, OOTResult::ootfracerr, OOTResult::ootfracsum, OOTResult::ootfracsumerr, and mathSSE::sqrt().

Referenced by ComputeOOTFractionvsFill(), and ComputeOOTFractionvsRun().

75  {
76 
77  if(perFill) {std::cout << "Processing fill " << run << std::endl;} else {std::cout << "Processing run " << run << std::endl;}
78 
79  char itpath[100];
80  char ootpath[100];
81  char etpath[100];
82  if(perFill) {
83  sprintf(itpath,"%s/fill_%d",itmodule,run);
84  sprintf(ootpath,"%s/fill_%d",ootmodule,run);
85  sprintf(etpath,"%s/fill_%d",etmodule,run);
86  }
87  else {
88  sprintf(itpath,"%s/run_%d",itmodule,run);
89  sprintf(ootpath,"%s/run_%d",ootmodule,run);
90  sprintf(etpath,"%s/run_%d",etmodule,run);
91  }
92 
93  OOTResult* res = new OOTResult;
94 
95  std::vector<int> filledbx = FillingSchemeFromProfile(ff,itpath,hname);
96  res->nfilledbx = filledbx.size();
97 
98  if(!perFill) {
99  std::vector<int> filledbxtest = FillingScheme(ff,etpath);
100  if(filledbx.size() != filledbxtest.size()) std::cout << "Inconsistency in number of filled BX "
101  << run << " " << filledbx.size() << " " << filledbxtest.size() << std::endl;
102  }
103 
104  TProfile* itmult = nullptr;
105  TProfile* ootmult = nullptr;
106  res->hratio = new TH1F("ratio","ratio",200,0.,2.);
107 
108  float rclzb=0;
109  float rclrandom=0;
110  float errclzb=0;
111  float errclrandom=0;
112  float nzb=0;
113  float nrandom=0;
114 
115  if(ff) {
116  if(ff->cd(itpath)) {
117  itmult = (TProfile*)gDirectory->Get(hname);
118  }
119  else { std::cout << "In time path is not ok" << std::endl;}
120  if(ff->cd(ootpath)) {
121  ootmult = (TProfile*)gDirectory->Get(hname);
122  }
123  else { std::cout << "out of time path is not ok" << std::endl;}
124  if(itmult && ootmult) {
125  ootmult->SetLineColor(kRed); ootmult->SetMarkerColor(kRed);
126  // ootmult->Draw();
127  // itmult->Draw("same");
128  for(std::vector<int>::const_iterator fbx=filledbx.begin();fbx!=filledbx.end();++fbx) {
129  nzb += itmult->GetBinEntries(*fbx);
130  nrandom += ootmult->GetBinEntries(*fbx+1);
131  }
132  for(std::vector<int>::const_iterator fbx=filledbx.begin();fbx!=filledbx.end();++fbx) {
133  if(nzb > 0 && nrandom > 0) {
134  rclzb += (itmult->GetBinContent(*fbx)*itmult->GetBinEntries(*fbx))/nzb;
135  errclzb += (itmult->GetBinError(*fbx)*itmult->GetBinEntries(*fbx))*(itmult->GetBinError(*fbx)*itmult->GetBinEntries(*fbx))/(nzb*nzb);
136  rclrandom += (ootmult->GetBinContent(*fbx+1)*ootmult->GetBinEntries(*fbx+1))/nrandom;
137  errclrandom += (ootmult->GetBinError(*fbx+1)*ootmult->GetBinEntries(*fbx+1))*(ootmult->GetBinError(*fbx+1)*ootmult->GetBinEntries(*fbx+1))
138  /(nrandom*nrandom);
139  }
140  if(itmult->GetBinContent(*fbx)==0) { std::cout << "No cluster in filled BX! " << *fbx << std::endl;}
141  else if(ootmult->GetBinEntries(*fbx+1)==0) {/* std::cout << "No entry in OOT BX " << *fbx+1 << std::endl; */}
142  else {
143  float rat= ootmult->GetBinContent(*fbx+1)/itmult->GetBinContent(*fbx);
144  res->hratio->Fill(rat);
145  }
146  }
147  }
148  else {std::cout<<"histograms not found"<<std::endl;}
149  }
150  else { std::cout << "Input file pointer is not ok" <<std::endl;}
151 
152  res->ngoodbx = res->hratio->GetEntries();
153 
154  if(nzb>0 && nrandom>0 &&rclzb>0 && rclrandom > 0) {
155  res->ootfracsum = rclrandom/rclzb;
156  res->ootfracsumerr = rclrandom/rclzb*sqrt(errclzb/(rclzb*rclzb)+errclrandom/(rclrandom*rclrandom));
157  }
158  if(res->ngoodbx) {
159  res->hratio->Fit("gaus","Q0","",.01,1.99);
160  if(res->hratio->GetFunction("gaus")) {
161  res->ootfrac = res->hratio->GetFunction("gaus")->GetParameter(1);
162  res->ootfracerr = res->hratio->GetFunction("gaus")->GetParError(1);
163  }
164  else {std::cout << "Missing fitting function" << std::endl;}
165  }
166  else {std::cout << "No filled BX or strange filling scheme" <<std::endl;}
167 
168  return res;
169 }
Definition: Electron.h:6
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< int > FillingSchemeFromProfile(TFile *ff, const char *path, const char *hname, const float thr)
std::vector< int > FillingScheme(TFile *ff, const char *path, const float thr)
OOTSummary* ComputeOOTFractionvsFill ( TFile *  ff,
const char *  itmodule,
const char *  ootmodule,
const char *  etmodule,
const char *  hname,
OOTSummary ootsumm = nullptr 
)

Definition at line 12 of file OOTMultiplicityPlotMacros.cc.

References stringResolutionProvider_cfi::bin, ComputeOOTFraction(), gather_cfg::cout, CommonAnalyzer::getFillList(), OOTSummary::hngoodbx, OOTSummary::hootfrac, OOTSummary::hootfracsum, OOTResult::hratio, mps_fire::i, OOTResult::ngoodbx, OOTResult::ootfrac, OOTResult::ootfracerr, OOTResult::ootfracsum, OOTResult::ootfracsumerr, writedatasetfile::runs, and jetUpdater_cfi::sort.

12  {
13 
14  if(ootsumm==nullptr) {
15  ootsumm = new OOTSummary;
16  }
17 
18  if(ff) {
19  CommonAnalyzer ca(ff,"",itmodule);
20  std::vector<unsigned int> runs = ca.getFillList();
21  std::sort(runs.begin(),runs.end());
22  for(unsigned int i=0;i<runs.size();++i) {
23 
24  char runlabel[100];
25  sprintf(runlabel,"%d",runs[i]);
26 
27  OOTResult* res = ComputeOOTFraction(ff,itmodule,ootmodule,etmodule,runs[i],hname,true);
28 
29  if(res->ngoodbx!=res->hratio->GetEntries()) std::cout << "Inconsistency in number of good bx" << std::endl;
30  ootsumm->hngoodbx->Fill(runlabel,res->ngoodbx);
31  int ibin = ootsumm->hootfracsum->Fill(runlabel,res->ootfracsum);
32  ootsumm->hootfracsum->SetBinError(ibin,res->ootfracsumerr);
33  int bin = ootsumm->hootfrac->Fill(runlabel,res->ootfrac);
34  ootsumm->hootfrac->SetBinError(bin,res->ootfracerr);
35  delete res;
36 
37  }
38  }
39  else {std::cout << "File is not ok" << std::endl; }
40 
41  return ootsumm;
42 }
Definition: Electron.h:6
OOTResult * ComputeOOTFraction(TFile *ff, const char *itmodule, const char *ootmodule, const char *etmodule, const int run, const char *hname, const bool &perFill)
bin
set the eta bin as selection string.
OOTSummary* ComputeOOTFractionvsRun ( TFile *  ff,
const char *  itmodule,
const char *  ootmodule,
const char *  etmodule,
const char *  hname,
OOTSummary ootsumm = nullptr 
)

Definition at line 43 of file OOTMultiplicityPlotMacros.cc.

References stringResolutionProvider_cfi::bin, ComputeOOTFraction(), gather_cfg::cout, CommonAnalyzer::getRunList(), OOTSummary::hngoodbx, OOTSummary::hootfrac, OOTSummary::hootfracsum, OOTResult::hratio, mps_fire::i, OOTResult::ngoodbx, OOTResult::ootfrac, OOTResult::ootfracerr, OOTResult::ootfracsum, OOTResult::ootfracsumerr, writedatasetfile::runs, and jetUpdater_cfi::sort.

43  {
44 
45  if(ootsumm==nullptr) {
46  ootsumm =new OOTSummary;
47  }
48 
49  if(ff) {
50  CommonAnalyzer ca(ff,"",itmodule);
51  std::vector<unsigned int> runs = ca.getRunList();
52  std::sort(runs.begin(),runs.end());
53  for(unsigned int i=0;i<runs.size();++i) {
54 
55  char runlabel[100];
56  sprintf(runlabel,"%d",runs[i]);
57 
58  OOTResult* res = ComputeOOTFraction(ff,itmodule,ootmodule,etmodule,runs[i],hname);
59 
60  if(res->ngoodbx!=res->hratio->GetEntries()) std::cout << "Inconsistency in number of good bx" << std::endl;
61  ootsumm->hngoodbx->Fill(runlabel,res->ngoodbx);
62  int ibin = ootsumm->hootfracsum->Fill(runlabel,res->ootfracsum);
63  ootsumm->hootfracsum->SetBinError(ibin,res->ootfracsumerr);
64  int bin = ootsumm->hootfrac->Fill(runlabel,res->ootfrac);
65  ootsumm->hootfrac->SetBinError(bin,res->ootfracerr);
66  delete res;
67 
68  }
69  }
70  else {std::cout << "File is not ok" << std::endl; }
71 
72  return ootsumm;
73 }
Definition: Electron.h:6
OOTResult * ComputeOOTFraction(TFile *ff, const char *itmodule, const char *ootmodule, const char *etmodule, const int run, const char *hname, const bool &perFill)
bin
set the eta bin as selection string.
std::vector<int> FillingScheme ( TFile *  ff,
const char *  path,
const float  thr = 0. 
)

Definition at line 171 of file OOTMultiplicityPlotMacros.cc.

References gather_cfg::cout, and mps_fire::i.

Referenced by ComputeOOTFraction().

171  {
172 
173  TH1F* bx=nullptr;
174  std::vector<int> filledbx;
175  if(ff) {
176  if(ff->cd(path)) {
177  bx = (TH1F*)gDirectory->Get("bx");
178  if(bx) {
179  // bx->Draw();
180  std::cout << "Number of entries " << bx->GetEntries() << " threshold " << thr/3564.*bx->GetEntries() << std::endl;
181  for(int i=1;i<bx->GetNbinsX()+1;++i) {
182  if(bx->GetBinContent(i)>thr/3564.*bx->GetEntries()) {
183  if(!filledbx.empty() && i == filledbx[filledbx.size()-1]+1) {
184  std::cout << "This is not a 50ns run ! " << std::endl;
185  filledbx.clear();
186  return filledbx;
187  }
188  filledbx.push_back(i);
189  }
190  }
191  }
192  else { std::cout << "Histogram not found" << std::endl;}
193  }
194  else { std::cout << "module path is not ok" << std::endl;}
195  }
196  else { std::cout << "Input file pointer is not ok" <<std::endl;}
197 
198  // std::cout << filledbx.size() << " filled bunch crossings" << std::endl;
199  // for(std::vector<int>::const_iterator fbx=filledbx.begin();fbx!=filledbx.end();++fbx) { std::cout << *fbx << std::endl;}
200  return filledbx;
201 }
std::vector<int> FillingSchemeFromProfile ( TFile *  ff,
const char *  path,
const char *  hname,
const float  thr = 0. 
)

Definition at line 202 of file OOTMultiplicityPlotMacros.cc.

References gather_cfg::cout, and mps_fire::i.

Referenced by ComputeOOTFraction().

202  {
203 
204  TProfile* bx=nullptr;
205  std::vector<int> filledbx;
206  if(ff) {
207  if(ff->cd(path)) {
208  bx = (TProfile*)gDirectory->Get(hname);
209  if(bx) {
210  // bx->Draw();
211  std::cout << "Number of entries " << bx->GetEntries() << " threshold " << thr/3564.*bx->GetEntries() << std::endl;
212  for(int i=1;i<bx->GetNbinsX()+1;++i) {
213  if(bx->GetBinEntries(i)>thr/3564.*bx->GetEntries()) {
214  if(!filledbx.empty() && i == filledbx[filledbx.size()-1]+1) {
215  std::cout << "This is not a 50ns run ! " << std::endl;
216  filledbx.clear();
217  return filledbx;
218  }
219  filledbx.push_back(i);
220  }
221  }
222  }
223  else { std::cout << "Histogram not found" << std::endl;}
224  }
225  else { std::cout << "module path is not ok" << std::endl;}
226  }
227  else { std::cout << "Input file pointer is not ok" <<std::endl;}
228 
229  // std::cout << filledbx.size() << " filled bunch crossings" << std::endl;
230  // for(std::vector<int>::const_iterator fbx=filledbx.begin();fbx!=filledbx.end();++fbx) { std::cout << *fbx << std::endl;}
231  return filledbx;
232 }