CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
CompareAlignments Class Reference

#include <CompareAlignments.h>

Public Member Functions

 CompareAlignments (const std::string &outPath)
 
void doComparison (TString namesandlabels, TString legendheader="", TString lefttitle="", TString righttitle="", PublicationStatus status=INTERNAL, bool bigtext=false)
 

Private Member Functions

void ColourStatsBoxes (TObjArray *hists)
 
void MergeRootfile (TDirectory *target, TList *sourcelist, TList *labellist, bool bigtext)
 
void nicePad (Int_t logx, Int_t logy)
 
void SetMinMaxRange (TObjArray *hists)
 

Private Attributes

TList * FileList
 
TList * LabelList
 
std::vector< std::string > lowestlevels
 
std::string outPath
 
std::vector< int > phases
 
TFile * Target
 
std::vector< int > theColors
 
std::vector< int > theStyles
 

Detailed Description

Definition at line 25 of file CompareAlignments.h.

Constructor & Destructor Documentation

◆ CompareAlignments()

CompareAlignments::CompareAlignments ( const std::string &  outPath)
inline

Definition at line 43 of file CompareAlignments.h.

43 : outPath(outPath) {}

Member Function Documentation

◆ ColourStatsBoxes()

void CompareAlignments::ColourStatsBoxes ( TObjArray *  hists)
private

Definition at line 301 of file CompareAlignments.cc.

References HLT_2023v12_cff::Class, h, compare::hists, dqmMemoryStats::stats, testProducerWithPsetDescEmpty_cfi::x1, testProducerWithPsetDescEmpty_cfi::x2, testProducerWithPsetDescEmpty_cfi::y1, and testProducerWithPsetDescEmpty_cfi::y2.

Referenced by MergeRootfile().

301  {
302  Double_t fStatsX1 = 0.85, fStatsX2 = 1., fStatsY1 = 0.77, fStatsY2 = 0.92;
303  // colours stats boxes like hists' line colors and moves the next to each other
304  if (!hists)
305  return;
306  Double_t x1 = fStatsX1, x2 = fStatsX2, y1 = fStatsY1, y2 = fStatsY2;
307  for (Int_t iH = 0; iH < hists->GetEntries(); ++iH) {
308  TH1 *h = static_cast<TH1 *>(hists->At(iH));
309  if (!h)
310  continue;
311  TObject *statObj = h->GetListOfFunctions()->FindObject("stats");
312  if (statObj && statObj->InheritsFrom(TPaveStats::Class())) {
313  TPaveStats *stats = static_cast<TPaveStats *>(statObj);
314  stats->SetLineColor(static_cast<TH1 *>(hists->At(iH))->GetLineColor());
315  stats->SetTextColor(static_cast<TH1 *>(hists->At(iH))->GetLineColor());
316  stats->SetFillColor(10);
317  stats->SetX1NDC(x1);
318  stats->SetX2NDC(x2);
319  stats->SetY1NDC(y1);
320  stats->SetY2NDC(y2);
321  y2 = y1 - 0.005; // shift down 2
322  y1 = y2 - (fStatsY2 - fStatsY1); // shift down 1
323  if (y1 < 0.) {
324  y1 = fStatsY1;
325  y2 = fStatsY2; // restart y-positions
326  x2 = x1 - 0.005; // shift left 2
327  x1 = x2 - (fStatsX2 - fStatsX1); // shift left 1
328  if (x1 < 0.) { // give up, start again:
329  x1 = fStatsX1, x2 = fStatsX2, y1 = fStatsY1, y2 = fStatsY2;
330  }
331  }
332  stats->DrawClone();
333  //} else if (gStyle->GetOptStat() != 0) { // failure in case changed in list via TExec....
334  //this->Warning("ColourStatsBoxes", "No stats found for %s", hists->At(iH)->GetName());
335  }
336  }
337 }
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4

◆ doComparison()

void CompareAlignments::doComparison ( TString  namesandlabels,
TString  legendheader = "",
TString  lefttitle = "",
TString  righttitle = "",
PublicationStatus  status = INTERNAL,
bool  bigtext = false 
)

Definition at line 3 of file CompareAlignments.cc.

References cms::cuda::assert(), gather_cfg::cout, FileList, mps_fire::i, LabelList, TkAlStyle::legendheader, lowestlevels, MergeRootfile(), NONE, outPath, phases, TkAlStyle::set(), mps_update::status, Target, theColors, and theStyles.

Referenced by merge().

8  {
9  std::cout << "Comparing using: >" << namesandlabels << "<" << std::endl;
10 
11  TkAlStyle::legendheader = legendheader;
12  TkAlStyle::set(status, NONE, lefttitle, righttitle);
13  gStyle->SetOptStat(111110);
14  gStyle->SetOptTitle(0);
15 
16  Target = TFile::Open((outPath + "/result.root").c_str(), "RECREATE");
17  FileList = new TList();
18  LabelList = new TList();
19 
20  int formatCounter = 1;
21  //TObjArray* stringarray = namesandlabels.Tokenize(",");
22  TObjArray *nameandlabelpairs = namesandlabels.Tokenize(",");
23  for (Int_t i = 0; i < nameandlabelpairs->GetEntries(); ++i) {
24  TObjArray *aFileLegPair = TString(nameandlabelpairs->At(i)->GetName()).Tokenize("=");
25 
26  if (aFileLegPair->GetEntries() == 2) {
27  TFile *currentFile = TFile::Open(aFileLegPair->At(0)->GetName());
28  if (currentFile != nullptr && !currentFile->IsZombie()) {
29  FileList->Add(currentFile); // 2
30  if (currentFile->Get("TrackerOfflineValidation/Pixel/P1PXBBarrel_1")) {
31  phases.push_back(1);
32  } else if (currentFile->Get("TrackerOfflineValidation/Pixel/TPBBarrel_1")) {
33  phases.push_back(0);
34  } else {
35  std::cout << "Unknown phase for file " << aFileLegPair->At(0)->GetName() << std::endl;
36  assert(false);
37  }
38  if (TString(aFileLegPair->At(1)->GetName()).Contains("|")) {
39  TObjArray *formatedLegendEntry = TString(aFileLegPair->At(1)->GetName()).Tokenize("|");
40  LabelList->Add(formatedLegendEntry->At(0));
41  if (formatedLegendEntry->GetEntries() > 1) {
42  theColors.push_back(atoi(formatedLegendEntry->At(1)->GetName()));
43 
44  if (formatedLegendEntry->GetEntries() > 2)
45  theStyles.push_back(atoi(formatedLegendEntry->At(2)->GetName()) % 100);
46  else
47  theStyles.push_back(formatCounter);
48  } else {
49  std::cout << "if you give a \"|\" in the legend name you will need to at least give a int for the color"
50  << std::endl;
51  }
52  formatCounter++;
53  } else {
54  LabelList->Add(aFileLegPair->At(1));
55  theColors.push_back(formatCounter);
56  theStyles.push_back(formatCounter);
57  formatCounter++;
58  }
59  } else {
60  std::cout << "Could not open: " << aFileLegPair->At(0)->GetName() << std::endl;
61  }
62  } else {
63  std::cout << "Please give file name and legend entry in the following form:\n"
64  << " filename1=legendentry1,filename2=legendentry2[|color[|style]]" << std::endl;
65  }
66  }
67 
68  // ************************************************************
69  // List of Files
70  //FileList->Add( TFile::Open("../test/AlignmentValidation_Elliptical.root") );
71  //FileList->Add( TFile::Open("../test/AlignmentValidation_10pb.root") ); // 2
72  //FileList->Add( TFile::Open("../test/AlignmentValidation_custom.root") ); // 2
73  // ************************************************************
74 
75  // put here the lowest level up to which you want to combine the
76  // histograms
77  lowestlevels.push_back("TPBLadder");
78  lowestlevels.push_back("TPEPanel");
79  lowestlevels.push_back("TIBHalfShell");
80  lowestlevels.push_back("TIDRing");
81  lowestlevels.push_back("TOBRod");
82  lowestlevels.push_back("TECSide");
83  // phase 1
84  // it checks each one independently so no harm in having
85  // both phase 0 and 1 together in the vector
86  lowestlevels.push_back("P1PXBLadder");
87  lowestlevels.push_back("P1PXECPanel");
88 
90 }
std::vector< int > theColors
assert(be >=bs)
static void set(const PublicationStatus status, const Era era=NONE, const TString customTitle="", const TString customRightTitle="")
Definition: TkAlStyle.cc:275
std::vector< int > theStyles
static TString legendheader
Definition: TkAlStyle.h:101
std::vector< std::string > lowestlevels
std::vector< int > phases
Definition: TkAlStyle.h:43
void MergeRootfile(TDirectory *target, TList *sourcelist, TList *labellist, bool bigtext)

◆ MergeRootfile()

void CompareAlignments::MergeRootfile ( TDirectory *  target,
TList *  sourcelist,
TList *  labellist,
bool  bigtext 
)
private

Definition at line 92 of file CompareAlignments.cc.

References HltBtagPostValidation_cff::c, DMR_cfg::cerr, HLT_2023v12_cff::Class, ColourStatsBoxes(), gather_cfg::cout, BTVHLTOfflineSource_cfi::dirname, TkAlStyle::drawStandardTitle(), createfilelist::int, submitPVResolutionJobs::key, TkAlStyle::legendheader, lowestlevels, SiStripPI::max, nicePad(), getGTfromDQMFile::obj, castor_dqm_sourceclient_file_cfg::path, phases, cms::alpakatools::detail::power(), submitPVResolutionJobs::q, pfClustersFromCombinedCaloHF_cfi::scale, SetMinMaxRange(), mps_update::status, AlCaHLTBitMon_QueryRunRegistry::string, filterCSVwithJSON::target, TkAlStyle::textSize, theColors, and theStyles.

Referenced by doComparison().

92  {
93  if (sourcelist->GetSize() == 0) {
94  std::cout << "Cowardly refuse to merge empty SourceList! " << std::endl;
95  return;
96  }
97 
98  TString path((char *)strstr(target->GetPath(), ":"));
99  path.Remove(0, 2);
100 
101  TFile *first_source = (TFile *)sourcelist->First();
102  TObjString *first_label = (TObjString *)labellist->First();
103 
104  first_source->cd(path);
105  TDirectory *current_sourcedir = gDirectory;
106  //gain time, do not add the objects in the list in memory
107  Bool_t status = TH1::AddDirectoryStatus();
108  TH1::AddDirectory(kFALSE);
109 
110  // loop over all keys in this directory
111 
112  TIter nextkey(current_sourcedir->GetListOfKeys());
113  TKey *key, *oldkey = nullptr;
114  while ((key = (TKey *)nextkey())) {
115  //keep only the highest cycle number for each key
116  if (oldkey && !strcmp(oldkey->GetName(), key->GetName()))
117  continue;
118 
119  // read object from first source file
120  first_source->cd(path);
121  TObject *obj = key->ReadObj();
122 
123  auto itphase = phases.begin();
124  int firstfilephase = *itphase;
125 
126  if (obj->IsA()->InheritsFrom(TH1::Class())) {
127  // descendant of TH1 -> merge it
128  TCanvas c(obj->GetName(), obj->GetName());
129  c.SetFillColor(10);
130 
131  bool is2d = false;
132  if (strstr(obj->ClassName(), "TH2") != nullptr)
133  is2d = true;
134  TH1 *h1 = static_cast<TH1 *>(obj);
135 
136  int q = 1;
137  TObjArray *histarray = new TObjArray;
138 
139  h1->SetLineStyle(theStyles.at(q - 1));
140  h1->SetLineWidth(2);
141 
142  h1->SetTitle("");
143 
144  h1->SetLineColor(theColors.at(q - 1));
145  h1->GetYaxis()->SetTitleOffset(1.5);
146  if (strstr(h1->GetName(), "summary") != nullptr)
147  h1->Draw("x0e1*H");
148  else if (is2d)
149  h1->Draw();
150  else
151  h1->Draw();
152 
153  double max = h1->GetMaximum();
154  double scale = 1;
155  if (max > 1000) {
156  int power = (int)(TMath::Log10(max / 100));
157  scale = 1 / TMath::Power(10, power);
158  h1->GetYaxis()->SetTitle(((TString(h1->GetYaxis()->GetTitle()) + " [#times 10^{") += (int)power) + "}]");
159  }
160  h1->Scale(scale);
161 
162  int nPlots = sourcelist->GetSize();
163  double legendY = 0.80;
164  if (nPlots > 3) {
165  legendY -= 0.01 * (nPlots - 3);
166  }
167  if (bigtext) {
168  legendY -= 0.05;
169  }
170  if (legendY < 0.6) {
171  std::cerr << "Warning: Huge legend!" << std::endl;
172  legendY = 0.6;
173  }
174  TLegend leg(0.17, legendY, 0.85, 0.88);
175  bool hasheader = (TkAlStyle::legendheader != "");
176  if (hasheader)
177  leg.SetHeader(TkAlStyle::legendheader);
178  if (bigtext)
179  leg.SetTextSize(TkAlStyle::textSize);
180  leg.AddEntry(h1, first_label->String().Data(), "L");
181  leg.SetBorderSize(0);
182  leg.SetFillColor(10);
183  // loop over all source files and add the content of the
184  // correspondant histogram to the one pointed to by "h1"
185  TFile *nextsource = (TFile *)sourcelist->After(first_source);
186  TObjString *nextlabel = (TObjString *)labellist->After(labellist->First());
187 
188  histarray->Add(h1);
189  while (nextsource) {
190  TKey *key2;
191  bool wrongphase = false;
192  ++itphase;
193 
194  if (firstfilephase != *itphase && path.Contains("TrackerOfflineValidation/Pixel")) {
195  //skip this one
196  key2 = nullptr;
197  wrongphase = true;
198  } else {
199  // make sure we are at the correct directory level by cd'ing to path
200  nextsource->cd(path);
201  key2 = (TKey *)gDirectory->GetListOfKeys()->FindObject(h1->GetName());
202  }
203  if (key2) {
204  ++q;
205  TH1 *h2 = (TH1 *)key2->ReadObj();
206 
207  if (!is2d) {
208  h2->SetLineStyle(theStyles.at(q - 1));
209  h2->SetLineWidth(2);
210  }
211 
212  h2->SetLineColor(theColors.at(q - 1));
213  std::stringstream newname;
214  newname << h2->GetName() << q;
215  h2->Scale(scale);
216 
217  h2->SetName(newname.str().c_str());
218  if (strstr(newname.str().c_str(), "summary") != nullptr)
219  h2->DrawClone("x0*He1sames");
220  else if (is2d)
221  h2->DrawClone("sames");
222  else
223  h2->DrawClone("sames");
224  leg.AddEntry(c.FindObject(h2->GetName()), nextlabel->String().Data(), "L");
225  histarray->Add(c.FindObject(h2->GetName()));
226  delete h2;
227  } else if (wrongphase) {
228  //nothing
229  } else {
230  std::cerr << "Histogram " << path << " is not present in file " << nextsource->GetName() << std::endl;
231  }
232 
233  nextsource = (TFile *)sourcelist->After(nextsource);
234  nextlabel = (TObjString *)labellist->After(nextlabel);
235  }
236  nicePad(0, 0);
237  leg.Draw();
239  c.Update();
240  if (strstr(h1->GetName(), "summary") == nullptr)
241  SetMinMaxRange(histarray);
242  ColourStatsBoxes(histarray);
243  target->cd();
244  c.Update();
245  c.Write();
246  histarray->Delete();
247 
248  } else if (obj->IsA()->InheritsFrom(TDirectory::Class())) {
249  // it's a subdirectory
250 
251  std::string dirname = obj->GetName();
252  for (std::vector<std::string>::const_iterator lowlevelit = lowestlevels.begin(),
253  lowlevelitend = lowestlevels.end();
254  lowlevelit != lowlevelitend;
255  ++lowlevelit)
256  if (dirname.find(*lowlevelit) != std::string::npos)
257  return;
258 
259  // create a new subdir of same name and title in the target file
260  target->cd();
261  TDirectory *newdir = target->mkdir(obj->GetName(), obj->GetTitle());
262 
263  // newdir is now the starting point of another round of merging
264  // newdir still knows its depth within the target file via
265  // GetPath(), so we can still figure out where we are in the recursion
266  MergeRootfile(newdir, sourcelist, labellist, bigtext);
267 
268  } else {
269  // object is of no type that we know or can handle
270  std::cout << "Unknown object type, name: " << obj->GetName() << " title: " << obj->GetTitle() << std::endl;
271  }
272 
273  } // while ( ( TKey *key = (TKey*)nextkey() ) )
274 
275  // save modifications to target file
276  target->SaveSelf(kTRUE);
277  TH1::AddDirectory(status);
278 }
void nicePad(Int_t logx, Int_t logy)
std::vector< int > theColors
void SetMinMaxRange(TObjArray *hists)
void ColourStatsBoxes(TObjArray *hists)
static double textSize
Definition: TkAlStyle.h:103
std::vector< int > theStyles
static void drawStandardTitle()
Definition: TkAlStyle.h:75
static TString legendheader
Definition: TkAlStyle.h:101
key
prepare the HTCondor submission files and eventually submit them
std::vector< std::string > lowestlevels
std::vector< int > phases
constexpr unsigned int power(unsigned int base, unsigned int exponent)
void MergeRootfile(TDirectory *target, TList *sourcelist, TList *labellist, bool bigtext)

◆ nicePad()

void CompareAlignments::nicePad ( Int_t  logx,
Int_t  logy 
)
private

Definition at line 280 of file CompareAlignments.cc.

Referenced by MergeRootfile().

280  {
281  /*
282  gPad->SetBottomMargin(0.10);
283  gPad->SetRightMargin(0.1);
284  gPad->SetLeftMargin(0.15);
285  gPad->SetTopMargin(0.15);
286  gPad->SetTickx(1);
287  gPad->SetTicky(1);
288 */
289  if (logy == 1) {
290  gPad->SetLogy();
291  } else {
292  gPad->SetLogy(0);
293  }
294  if (logx == 1) {
295  gPad->SetLogx();
296  } else {
297  gPad->SetLogx(0);
298  }
299 }

◆ SetMinMaxRange()

void CompareAlignments::SetMinMaxRange ( TObjArray *  hists)
private

Definition at line 339 of file CompareAlignments.cc.

References h, compare::hists, mps_fire::i, SiStripPI::max, and SiStripPI::min.

Referenced by MergeRootfile().

339  {
340  Double_t min = 100000;
341  Double_t max = -100000;
342  for (Int_t iH = 0; iH < hists->GetEntries(); ++iH) {
343  TH1 *h = static_cast<TH1 *>(hists->At(iH));
344  if (!h)
345  continue;
346  for (int i = 1; i <= h->GetNbinsX(); ++i) {
347  if (h->GetBinContent(i) + h->GetBinError(i) > max)
348  max = h->GetBinContent(i) + h->GetBinError(i);
349  if (h->GetBinContent(i) - h->GetBinError(i) < min)
350  min = h->GetBinContent(i) - h->GetBinError(i);
351  }
352  }
353 
354  TH1 *h_first = static_cast<TH1 *>(hists->At(0));
355  h_first->SetMaximum(max * 1.3);
356  if (min == 0.) {
357  min = -1111;
358  h_first->SetMinimum(min);
359  } else {
360  h_first->SetMinimum(min - min * 0.1);
361  }
362 }
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4

Member Data Documentation

◆ FileList

TList* CompareAlignments::FileList
private

Definition at line 27 of file CompareAlignments.h.

Referenced by doComparison().

◆ LabelList

TList* CompareAlignments::LabelList
private

Definition at line 28 of file CompareAlignments.h.

Referenced by doComparison().

◆ lowestlevels

std::vector<std::string> CompareAlignments::lowestlevels
private

Definition at line 30 of file CompareAlignments.h.

Referenced by doComparison(), and MergeRootfile().

◆ outPath

std::string CompareAlignments::outPath
private

Definition at line 35 of file CompareAlignments.h.

Referenced by doComparison().

◆ phases

std::vector<int> CompareAlignments::phases
private

Definition at line 33 of file CompareAlignments.h.

Referenced by doComparison(), and MergeRootfile().

◆ Target

TFile* CompareAlignments::Target
private

Definition at line 29 of file CompareAlignments.h.

Referenced by doComparison().

◆ theColors

std::vector<int> CompareAlignments::theColors
private

Definition at line 31 of file CompareAlignments.h.

Referenced by doComparison(), and MergeRootfile().

◆ theStyles

std::vector<int> CompareAlignments::theStyles
private

Definition at line 32 of file CompareAlignments.h.

Referenced by doComparison(), and MergeRootfile().