CMS 3D CMS Logo

moduleOccupancyPlots.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <fstream>
3 #include <sstream>
4 #include <string>
5 
6 #include "TFile.h"
7 #include "TH1D.h"
8 #include "TDirectoryFile.h"
9 #include "TCanvas.h"
10 
11 #define NOISEPREFIX "Profile_NoiseFromCondDB__det__"
12 #define PEDESTALPREFIX "Profile_PedestalFromCondDB__det__"
13 #define OCCUPANCYPREFIX "ClusterDigiPosition__det__"
14 
15 void printPlot(TH1D* hist, char* prefix, char* postfix);
16 int getChannelNumber(int ival);
17 
18 int main(int argc, char* argv[]) {
19  char* rootfilename;
20  char* modulelistname;
21  int pnbits;
22  char* prefix;
23  char* postfix;
24 
25  bool debugPrint = false;
26  bool imagePrint = false;
27 
28  if (argc == 6) {
29  rootfilename = argv[1];
30  modulelistname = argv[2];
31  pnbits = atoi(argv[3]);
32  prefix = argv[4];
33  postfix = argv[5];
34  } else {
35  std::cout << "Wrong number of parameters " << argc << std::endl;
36  return 1;
37  }
38 
39  if (debugPrint)
40  std::cout << "ready to go " << rootfilename << ' ' << modulelistname << std::endl;
41 
42  TFile* rootfile = TFile::Open(rootfilename, "READ");
43  if (!rootfile) {
44  std::cout << "Problems with input root file" << std::endl;
45  return 2;
46  }
47  int detid;
48  std::ifstream modulelist(modulelistname);
49 
50  std::stringstream outrootfilename;
51  outrootfilename << prefix << "SummaryFile" << postfix << ".root";
52  TFile* outrootfile = TFile::Open(outrootfilename.str().c_str(), "RECREATE");
53 
54  TH1D* th_summary = nullptr;
55  Double_t TotalEvents = 0.0;
56  Double_t TotalDigis = 0.0;
57 
58  TH1D* TotEvents = new TH1D("TotEvents", "TotEvents", 1, 0, 1);
59 
60  if (pnbits & 4) {
61  TDirectoryFile* tdir = (TDirectoryFile*)rootfile->FindObjectAny("AlCaReco");
62  TH1D* hist = (TH1D*)tdir->FindObjectAny("TotalNumberOfCluster__TIB");
63  if (hist) {
64  TotalEvents = hist->GetEntries();
65  TotEvents->SetBinContent(1, TotalEvents);
66  TotEvents->Write();
67  }
68  }
69  while (modulelist >> detid) {
70  if (debugPrint)
71  std::cout << " ready to go with detid " << detid << " " << pnbits << std::endl;
72  // bit 0: noise
73  if (pnbits & 1) {
74  std::stringstream histoname;
75  histoname << NOISEPREFIX << detid;
76  if (debugPrint)
77  std::cout << " ready to go with histogram " << histoname.str() << std::endl;
78  TH1D* hist = (TH1D*)rootfile->FindObjectAny(histoname.str().c_str());
79  if (hist) {
80  if (debugPrint)
81  std::cout << histoname.str() << " found!" << std::endl;
82  if (imagePrint)
84  hist->Write();
85  } else {
86  std::cout << histoname.str() << " NOT found..." << std::endl;
87  }
88  }
89  // bit 1: pedestal
90  if (pnbits & 2) {
91  std::stringstream histoname;
92  histoname << PEDESTALPREFIX << detid;
93  if (debugPrint)
94  std::cout << " ready to go with histogram " << histoname.str() << std::endl;
95  TH1D* hist = (TH1D*)rootfile->FindObjectAny(histoname.str().c_str());
96  if (hist) {
97  if (debugPrint)
98  std::cout << histoname.str() << " found!" << std::endl;
99  if (imagePrint)
101  hist->Write();
102  } else {
103  std::cout << histoname.str() << " NOT found..." << std::endl;
104  }
105  }
106  // bit 2: Occupancy
107  if (pnbits & 4) {
108  std::stringstream histoname;
109  histoname << OCCUPANCYPREFIX << detid;
110  std::string SummaryName = "ClusterDigiPosition__det__Summary";
111 
112  if (debugPrint)
113  std::cout << " ready to go with histogram " << histoname.str() << std::endl;
114  if (th_summary == nullptr)
115  th_summary = new TH1D(SummaryName.c_str(), SummaryName.c_str(), 768, 0.5, 768.5);
116  TH1D* hist = (TH1D*)rootfile->FindObjectAny(histoname.str().c_str());
117  if (hist) {
118  if (debugPrint)
119  std::cout << histoname.str() << " found!" << hist->GetEntries() << std::endl;
120  for (int i = 1; i < hist->GetNbinsX() + 1; i++) {
121  th_summary->AddBinContent(i, hist->GetBinContent(i));
122  TotalDigis += hist->GetBinContent(i);
123  }
124  hist->SetLineColor(2);
125  if (imagePrint)
127  hist->Write();
128  } else {
129  std::cout << histoname.str() << " NOT found..." << std::endl;
130  }
131  }
132  // bit 3 : reorder
133  if (pnbits & 8) {
134  std::stringstream histoname;
135  histoname << OCCUPANCYPREFIX << detid;
136  if (debugPrint)
137  std::cout << " ready to go with histogram " << histoname.str() << std::endl;
138  TH1D* hist = (TH1D*)rootfile->FindObjectAny(histoname.str().c_str());
139  if (hist) {
140  if (debugPrint)
141  std::cout << histoname.str() << " found!" << std::endl;
142  std::stringstream histoname_reorder;
143  histoname_reorder << OCCUPANCYPREFIX << "_reorder_" << detid;
144  TH1D* hist_reorder = new TH1D(histoname_reorder.str().c_str(),
145  histoname_reorder.str().c_str(),
146  hist->GetXaxis()->GetNbins(),
147  hist->GetXaxis()->GetXmin(),
148  hist->GetXaxis()->GetXmax());
149  for (int i = 0; i < hist_reorder->GetNbinsX(); i++) {
150  int chan = getChannelNumber(i);
151  hist_reorder->SetBinContent(i + 1, hist->GetBinContent(chan));
152  }
153  hist->Write();
154  hist_reorder->Write();
155  } else {
156  std::cout << histoname.str() << " NOT found..." << std::endl;
157  }
158  }
159  }
160  if (th_summary) {
161  std::string fname = rootfilename;
162  std::string run = fname.substr(fname.find("R000") + 4, 6);
163  if (TotalEvents) {
164  Double_t fac = 1.0 / TotalEvents;
165  th_summary->Scale(fac);
166  std::cout << " Run Number " << run << " Events " << TotalEvents << " Total # of Digis " << TotalDigis
167  << " Av. Digis " << TotalDigis * fac << std::endl;
168  th_summary->SetEntries(TotalDigis * fac);
169  }
170  th_summary->Write();
171  printPlot(th_summary, prefix, postfix);
172  }
173 
174  outrootfile->Close();
175 
176  return 0;
177 }
178 
179 void printPlot(TH1D* hist, char* prefix, char* postfix) {
180  TCanvas* cc = new TCanvas;
181  hist->Draw();
182  std::stringstream filename;
183  filename << prefix << hist->GetName() << postfix << ".png";
184  cc->Print(filename.str().c_str());
185  delete cc;
186 }
187 int getChannelNumber(int ival) {
188  int chan = int(32 * fmod(int(fmod(ival, 256.) / 2.), 4.) + 8 * int(int(fmod(ival, 256.) / 2.) / 4.) -
189  31 * int(int(fmod(ival, 256.) / 2.) / 16.) + fmod(fmod(ival, 256.), 2.) * 128 + int(ival / 256) * 256);
190  return chan;
191 }
void printPlot(TH1D *hist, char *prefix, char *postfix)
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
int getChannelNumber(int ival)
#define NOISEPREFIX
#define OCCUPANCYPREFIX
chan
lumi = TPaveText(lowX+0.38, lowY+0.061, lowX+0.45, lowY+0.161, "NDC") lumi.SetBorderSize( 0 ) lumi...
string fname
main script
#define PEDESTALPREFIX
int main(int argc, char *argv[])